[10087] | 1 | ! |
---|
| 2 | ! $Id: modupdate.F 779 2007-12-22 17:04:17Z rblod $ |
---|
| 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 | !> Module Agrif_User_Interpolation |
---|
| 25 | !! |
---|
| 26 | !! This module contains procedures used by AGRIF users to do an interpolation |
---|
| 27 | !! either within the domain or at boundary or on the whole domain. |
---|
| 28 | !--------------------------------------------------------------------------------------------------- |
---|
| 29 | ! |
---|
| 30 | module Agrif_User_Interpolation |
---|
| 31 | ! |
---|
| 32 | use Agrif_Boundary |
---|
| 33 | use Agrif_user_variables |
---|
| 34 | ! |
---|
| 35 | Implicit none |
---|
| 36 | |
---|
| 37 | contains |
---|
| 38 | ! |
---|
| 39 | !=================================================================================================== |
---|
| 40 | ! subroutine Agrif_Set_bc |
---|
| 41 | !> This subroutine is used to specify the variable we want to interpolate and the boundary location of |
---|
| 42 | !! the current grid we want to interpolate (boundary interpolation). |
---|
| 43 | !--------------------------------------------------------------------------------------------------- |
---|
| 44 | subroutine Agrif_Set_bc ( tabvarsindic, bcinfsup, Interpolationshouldbemade ) |
---|
| 45 | !--------------------------------------------------------------------------------------------------- |
---|
| 46 | integer, intent(in) :: tabvarsindic !< index of the variable in tabvars |
---|
| 47 | integer, dimension(2), intent(in) :: bcinfsup !< location |
---|
| 48 | logical, optional, intent(in) :: Interpolationshouldbemade !< interpolation should be made |
---|
| 49 | ! |
---|
| 50 | integer :: indic ! indice of the variable in tabvars |
---|
| 51 | type(Agrif_Variable), pointer :: var |
---|
| 52 | ! |
---|
| 53 | var => Agrif_Search_Variable(Agrif_Curgrid,tabvarsindic) |
---|
| 54 | if (.not.associated(var)) return ! Grand mother grid case |
---|
| 55 | ! |
---|
| 56 | if ( Agrif_Curgrid % fixedrank /= 0 ) then |
---|
| 57 | if ( .not.associated(var % oldvalues2D) ) then |
---|
| 58 | allocate(var % oldvalues2D(2,1)) |
---|
| 59 | var % interpIndex = -1 |
---|
| 60 | var % oldvalues2D = 0. |
---|
| 61 | endif |
---|
| 62 | if ( present(Interpolationshouldbemade) ) then |
---|
| 63 | var % Interpolationshouldbemade = Interpolationshouldbemade |
---|
| 64 | endif |
---|
| 65 | endif |
---|
| 66 | ! |
---|
| 67 | var % bcinf = bcinfsup(1) |
---|
| 68 | var % bcsup = bcinfsup(2) |
---|
| 69 | !--------------------------------------------------------------------------------------------------- |
---|
| 70 | end subroutine Agrif_Set_bc |
---|
| 71 | !=================================================================================================== |
---|
| 72 | ! |
---|
| 73 | !=================================================================================================== |
---|
| 74 | ! subroutine Agrif_Set_interp |
---|
| 75 | !> This procedure specify the type of interpolation in case of interpolation within the domain. |
---|
| 76 | !! The index of the variable we want to interpolate is a profile defined with the subroutine Agrif_Declare_Variable. |
---|
| 77 | ! |
---|
| 78 | ! Interp1 indicates the type of interpolation in the first dimension |
---|
| 79 | !--------------------------------------------------------------------------------------------------- |
---|
| 80 | subroutine Agrif_Set_interp ( tabvarsindic, interp, interp1, interp2, interp3 , interp4) |
---|
| 81 | !--------------------------------------------------------------------------------------------------- |
---|
| 82 | integer, intent(in) :: tabvarsindic !< index of the variable in tabvars |
---|
| 83 | integer, optional, intent(in) :: interp, interp1, interp2, interp3, interp4 !< The type of interpolation |
---|
| 84 | ! |
---|
| 85 | integer :: indic ! index of the variable in tabvars |
---|
| 86 | type(Agrif_Variable), pointer :: var |
---|
| 87 | ! |
---|
| 88 | var => Agrif_Search_Variable(Agrif_Curgrid,tabvarsindic) |
---|
| 89 | if (.not.associated(var)) return ! Grand mother grid case |
---|
| 90 | ! |
---|
| 91 | var % type_interp = Agrif_Constant |
---|
| 92 | ! |
---|
| 93 | if (present(interp)) var % type_interp = interp |
---|
| 94 | if (present(interp1)) var % type_interp(1) = interp1 |
---|
| 95 | if (present(interp2)) var % type_interp(2) = interp2 |
---|
| 96 | if (present(interp3)) var % type_interp(3) = interp3 |
---|
| 97 | if (present(interp4)) var % type_interp(4) = interp4 |
---|
| 98 | !--------------------------------------------------------------------------------------------------- |
---|
| 99 | end subroutine Agrif_Set_interp |
---|
| 100 | !=================================================================================================== |
---|
| 101 | ! |
---|
| 102 | !=================================================================================================== |
---|
| 103 | ! subroutine Agrif_Set_bcinterp |
---|
| 104 | !> This subroutine is used to indicate the type of interpolation in case of boundary interpolation. |
---|
| 105 | !! The index of the variable we want to interpolate is a profile defined with the subroutine Agrif_Declare_Variable. |
---|
| 106 | ! |
---|
| 107 | !> Interp1 indicates the type of interpolation in the first dimension |
---|
| 108 | !! Interp21 indicates the first variable of the second dimension |
---|
| 109 | !--------------------------------------------------------------------------------------------------- |
---|
| 110 | subroutine Agrif_Set_bcinterp ( tabvarsindic, interp, interp1, interp2, interp3, interp4, & |
---|
| 111 | interp11, interp12, interp21, interp22 ) |
---|
| 112 | !--------------------------------------------------------------------------------------------------- |
---|
| 113 | INTEGER, intent(in) :: tabvarsindic !< index of the variable in tabvars |
---|
| 114 | INTEGER, OPTIONAL, intent(in) :: interp, interp1, interp2, interp3, interp4 !< type of interpolation |
---|
| 115 | INTEGER, OPTIONAL, intent(in) :: interp11, interp12, interp21, interp22 !< dimension we want to interpolate |
---|
| 116 | ! |
---|
| 117 | INTEGER :: indic ! index of the variable in tabvars |
---|
| 118 | TYPE(Agrif_Variable), pointer :: var |
---|
| 119 | ! |
---|
| 120 | var => Agrif_Search_Variable(Agrif_Curgrid,tabvarsindic) |
---|
| 121 | ! |
---|
| 122 | var % type_interp_bc = Agrif_Constant |
---|
| 123 | ! |
---|
| 124 | if (present(interp)) var % type_interp_bc = interp |
---|
| 125 | if (present(interp1)) var % type_interp_bc(:,1) = interp1 |
---|
| 126 | if (present(interp11)) var % type_interp_bc(1,1) = interp11 |
---|
| 127 | if (present(interp12)) var % type_interp_bc(1,2) = interp12 |
---|
| 128 | if (present(interp2)) var % type_interp_bc(:,2) = interp2 |
---|
| 129 | if (present(interp21)) var % type_interp_bc(2,1) = interp21 |
---|
| 130 | if (present(interp22)) var % type_interp_bc(2,2) = interp22 |
---|
| 131 | if (present(interp3)) var % type_interp_bc(:,3) = interp3 |
---|
| 132 | if (present(interp4)) var % type_interp_bc(:,4) = interp4 |
---|
| 133 | !--------------------------------------------------------------------------------------------------- |
---|
| 134 | end subroutine Agrif_Set_bcinterp |
---|
| 135 | !=================================================================================================== |
---|
| 136 | ! |
---|
| 137 | !=================================================================================================== |
---|
| 138 | ! subroutine Agrif_Init_variable |
---|
| 139 | !> This subroutine is used for interpolation over the whole domain (within the domain and boundaries). |
---|
| 140 | !! It specifies the variable we want to interpolate and define through an input procedure (procname) |
---|
| 141 | !! what we do either on parent or child grid. The index of the variable we want to interpolate is a profile |
---|
| 142 | !! defined with the subroutine Agrif_Declare_Variable. |
---|
| 143 | !--------------------------------------------------------------------------------------------------- |
---|
| 144 | subroutine Agrif_Init_variable ( tabvarsindic, procname ) |
---|
| 145 | !--------------------------------------------------------------------------------------------------- |
---|
| 146 | INTEGER, intent(in) :: tabvarsindic !< index of the variable in tabvars |
---|
| 147 | procedure() :: procname !< Data recovery procedure written by users |
---|
| 148 | ! |
---|
| 149 | if ( Agrif_Curgrid%level <= 0 ) return |
---|
| 150 | ! |
---|
| 151 | call Agrif_Interp_variable(tabvarsindic, procname) |
---|
| 152 | call Agrif_Bc_variable(tabvarsindic, procname, 1.) |
---|
| 153 | !--------------------------------------------------------------------------------------------------- |
---|
| 154 | end subroutine Agrif_Init_variable |
---|
| 155 | !=================================================================================================== |
---|
| 156 | ! |
---|
| 157 | !=================================================================================================== |
---|
| 158 | ! subroutine Agrif_Bc_variable |
---|
| 159 | !> This subroutine is used for boundary interpolation. It specifies the variable we want to interpolate |
---|
| 160 | !! and define through an input procedure (procname) what we do either on parent or child grid. |
---|
| 161 | !! The index of the variable we want to interpolate is a profile defined with the subroutine Agrif_Declare_Variable. |
---|
| 162 | !--------------------------------------------------------------------------------------------------- |
---|
| 163 | subroutine Agrif_Bc_variable ( tabvarsindic, procname, calledweight ) |
---|
| 164 | !--------------------------------------------------------------------------------------------------- |
---|
| 165 | integer, intent(in) :: tabvarsindic !< index of the variable in tabvars |
---|
| 166 | procedure() :: procname !< Data recovery procedure written by users |
---|
| 167 | real, optional, intent(in) :: calledweight !< weight of interpolation |
---|
| 168 | ! |
---|
| 169 | real :: weight |
---|
| 170 | logical :: pweight |
---|
| 171 | integer :: indic |
---|
| 172 | integer :: nbdim |
---|
| 173 | type(Agrif_Variable), pointer :: root_var |
---|
| 174 | type(Agrif_Variable), pointer :: parent_var |
---|
| 175 | type(Agrif_Variable), pointer :: child_var |
---|
| 176 | type(Agrif_Variable), pointer :: child_tmp ! Temporary variable on the child grid |
---|
| 177 | integer :: i |
---|
| 178 | integer,dimension(7) :: lb, ub |
---|
| 179 | ! |
---|
| 180 | if ( Agrif_Curgrid%level <= 0 ) return |
---|
| 181 | ! |
---|
| 182 | ! |
---|
| 183 | if ( present(calledweight) ) then |
---|
| 184 | weight = calledweight |
---|
| 185 | pweight = .true. |
---|
| 186 | else |
---|
| 187 | weight = 0. |
---|
| 188 | pweight = .false. |
---|
| 189 | endif |
---|
| 190 | ! |
---|
| 191 | child_var => Agrif_Search_Variable(Agrif_Curgrid,tabvarsindic) |
---|
| 192 | parent_var => child_var % parent_var |
---|
| 193 | root_var => child_var % root_var |
---|
| 194 | ! |
---|
| 195 | nbdim = root_var % nbdim |
---|
| 196 | ! |
---|
| 197 | do i=1,nbdim |
---|
| 198 | if (root_var%coords(i) == 0) then |
---|
| 199 | lb(i) = parent_var%lb(i) |
---|
| 200 | ub(i) = parent_var%ub(i) |
---|
| 201 | else |
---|
| 202 | lb(i) = child_var%lb(i) |
---|
| 203 | ub(i) = child_var%ub(i) |
---|
| 204 | endif |
---|
| 205 | enddo |
---|
| 206 | |
---|
| 207 | select case( nbdim ) |
---|
| 208 | case(1) |
---|
| 209 | allocate(parray1(lb(1):ub(1))) |
---|
| 210 | case(2) |
---|
| 211 | allocate(parray2(lb(1):ub(1), & |
---|
| 212 | lb(2):ub(2) )) |
---|
| 213 | case(3) |
---|
| 214 | allocate(parray3(lb(1):ub(1), & |
---|
| 215 | lb(2):ub(2), & |
---|
| 216 | lb(3):ub(3) )) |
---|
| 217 | case(4) |
---|
| 218 | allocate(parray4(lb(1):ub(1), & |
---|
| 219 | lb(2):ub(2), & |
---|
| 220 | lb(3):ub(3), & |
---|
| 221 | lb(4):ub(4) )) |
---|
| 222 | case(5) |
---|
| 223 | allocate(parray5(lb(1):ub(1), & |
---|
| 224 | lb(2):ub(2), & |
---|
| 225 | lb(3):ub(3), & |
---|
| 226 | lb(4):ub(4), & |
---|
| 227 | lb(5):ub(5) )) |
---|
| 228 | case(6) |
---|
| 229 | allocate(parray6(lb(1):ub(1), & |
---|
| 230 | lb(2):ub(2), & |
---|
| 231 | lb(3):ub(3), & |
---|
| 232 | lb(4):ub(4), & |
---|
| 233 | lb(5):ub(5), & |
---|
| 234 | lb(6):ub(6) )) |
---|
| 235 | end select |
---|
| 236 | ! |
---|
| 237 | ! Create temporary child variable |
---|
| 238 | allocate(child_tmp) |
---|
| 239 | ! |
---|
| 240 | child_tmp % root_var => root_var |
---|
| 241 | child_tmp % oldvalues2D => child_var % oldvalues2D |
---|
| 242 | ! |
---|
| 243 | ! Index indicating if a space interpolation is necessary |
---|
| 244 | child_tmp % interpIndex = child_var % interpIndex |
---|
| 245 | child_tmp % list_interp => child_var % list_interp |
---|
| 246 | child_tmp % Interpolationshouldbemade = child_var % Interpolationshouldbemade |
---|
| 247 | ! |
---|
| 248 | child_tmp % point = child_var % point |
---|
| 249 | child_tmp % lb = child_var % lb |
---|
| 250 | child_tmp % ub = child_var % ub |
---|
| 251 | |
---|
| 252 | child_tmp % lubglob = child_var % lubglob |
---|
| 253 | |
---|
| 254 | ! |
---|
| 255 | child_tmp % bcinf = child_var % bcinf |
---|
| 256 | child_tmp % bcsup = child_var % bcsup |
---|
| 257 | ! |
---|
| 258 | child_tmp % childarray = child_var % childarray |
---|
| 259 | child_tmp % memberin = child_var % memberin |
---|
| 260 | ! |
---|
| 261 | call Agrif_CorrectVariable(parent_var, child_tmp, pweight, weight, procname) |
---|
| 262 | ! |
---|
| 263 | child_var % childarray = child_tmp % childarray |
---|
| 264 | child_var % memberin = child_tmp % memberin |
---|
| 265 | ! |
---|
| 266 | child_var % oldvalues2D => child_tmp % oldvalues2D |
---|
| 267 | child_var % list_interp => child_tmp % list_interp |
---|
| 268 | ! |
---|
| 269 | child_var % interpIndex = child_tmp % interpIndex |
---|
| 270 | ! |
---|
| 271 | deallocate(child_tmp) |
---|
| 272 | ! |
---|
| 273 | select case( nbdim ) |
---|
| 274 | case(1); deallocate(parray1) |
---|
| 275 | case(2); deallocate(parray2) |
---|
| 276 | case(3); deallocate(parray3) |
---|
| 277 | case(4); deallocate(parray4) |
---|
| 278 | case(5); deallocate(parray5) |
---|
| 279 | case(6); deallocate(parray6) |
---|
| 280 | end select |
---|
| 281 | !--------------------------------------------------------------------------------------------------- |
---|
| 282 | end subroutine Agrif_Bc_variable |
---|
| 283 | !=================================================================================================== |
---|
| 284 | ! |
---|
| 285 | !=================================================================================================== |
---|
| 286 | ! subroutine Agrif_Interp_variable |
---|
| 287 | !> This subroutine is used for interpolation within the domain. It specifies the variable we want to interpolate |
---|
| 288 | !! and define through an input procedure (procname) what we do either on parent or child grid. |
---|
| 289 | !! The index of the variable we want to interpolate is a profile defined with the subroutine Agrif_Declare_Variable. |
---|
| 290 | !--------------------------------------------------------------------------------------------------- |
---|
| 291 | subroutine Agrif_Interp_variable ( tabvarsindic, procname ) |
---|
| 292 | !--------------------------------------------------------------------------------------------------- |
---|
| 293 | integer, intent(in) :: tabvarsindic !< index of the variable in tabvars |
---|
| 294 | procedure() :: procname !< Data recovery procedure written by users |
---|
| 295 | ! |
---|
| 296 | integer :: nbdim |
---|
| 297 | integer :: indic ! indice of the variable in tabvars |
---|
| 298 | logical :: torestore |
---|
| 299 | type(Agrif_Variable), pointer :: root_var |
---|
| 300 | type(Agrif_Variable), pointer :: parent_var ! Variable on the parent grid |
---|
| 301 | type(Agrif_Variable), pointer :: child_var ! Variable on the parent grid |
---|
| 302 | type(Agrif_Variable), pointer :: child_tmp ! Temporary variable on the child grid |
---|
| 303 | ! |
---|
| 304 | if ( Agrif_Curgrid%level <= 0 ) return |
---|
| 305 | ! |
---|
| 306 | |
---|
| 307 | child_var => Agrif_Search_Variable(Agrif_Curgrid,tabvarsindic) |
---|
| 308 | parent_var => child_var % parent_var |
---|
| 309 | root_var => child_var % root_var |
---|
| 310 | |
---|
| 311 | ! |
---|
| 312 | nbdim = root_var % nbdim |
---|
| 313 | torestore = root_var % restore |
---|
| 314 | ! |
---|
| 315 | allocate(child_tmp) |
---|
| 316 | ! |
---|
| 317 | child_tmp % root_var => root_var |
---|
| 318 | child_tmp % nbdim = root_var % nbdim |
---|
| 319 | child_tmp % point = child_var % point |
---|
| 320 | child_tmp % lb = child_var % lb |
---|
| 321 | child_tmp % ub = child_var % ub |
---|
| 322 | child_tmp % lubglob = child_var % lubglob |
---|
| 323 | child_tmp % interpIndex = child_var % interpIndex |
---|
| 324 | child_tmp % list_interp => child_var % list_interp |
---|
| 325 | child_tmp % Interpolationshouldbemade = child_var % Interpolationshouldbemade |
---|
| 326 | ! |
---|
| 327 | if ( torestore ) then |
---|
| 328 | select case( nbdim ) |
---|
| 329 | case(1) |
---|
| 330 | parray1 = child_var % array1 |
---|
| 331 | child_tmp % restore1D => child_var % restore1D |
---|
| 332 | case(2) |
---|
| 333 | parray2 = child_var % array2 |
---|
| 334 | child_tmp % restore2D => child_var % restore2D |
---|
| 335 | case(3) |
---|
| 336 | parray3 = child_var % array3 |
---|
| 337 | child_tmp % restore3D => child_var % restore3D |
---|
| 338 | case(4) |
---|
| 339 | parray4 = child_var % array4 |
---|
| 340 | child_tmp % restore4D => child_var % restore4D |
---|
| 341 | case(5) |
---|
| 342 | parray5 = child_var % array5 |
---|
| 343 | child_tmp % restore5D => child_var % restore5D |
---|
| 344 | case(6) |
---|
| 345 | parray6 = child_var % array6 |
---|
| 346 | child_tmp % restore6D => child_var % restore6D |
---|
| 347 | end select |
---|
| 348 | endif |
---|
| 349 | ! |
---|
| 350 | call Agrif_InterpVariable(parent_var, child_tmp, torestore, procname) |
---|
| 351 | ! |
---|
| 352 | child_var % list_interp => child_tmp % list_interp |
---|
| 353 | ! |
---|
| 354 | deallocate(child_tmp) |
---|
| 355 | !--------------------------------------------------------------------------------------------------- |
---|
| 356 | end subroutine Agrif_Interp_variable |
---|
| 357 | !=================================================================================================== |
---|
| 358 | ! |
---|
| 359 | |
---|
| 360 | end module Agrif_User_Interpolation |
---|