[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_Variables |
---|
| 25 | !! |
---|
| 26 | !! This module contains a procedure which will help to define a new variable which will |
---|
| 27 | !! be used for interpolation/update. |
---|
| 28 | !---------------------------------------------------------------------------------------------------- |
---|
| 29 | ! |
---|
| 30 | module Agrif_User_Variables |
---|
| 31 | ! |
---|
| 32 | use Agrif_User_Functions |
---|
| 33 | ! |
---|
| 34 | implicit none |
---|
| 35 | ! |
---|
| 36 | contains |
---|
| 37 | ! |
---|
| 38 | !=================================================================================================== |
---|
| 39 | ! subroutine Agrif_Declare_Variable |
---|
| 40 | ! |
---|
| 41 | !> This subroutine is used to declare a profile (new variable) which will be used for interpolation/update |
---|
| 42 | !--------------------------------------------------------------------------------------------------- |
---|
| 43 | subroutine Agrif_Declare_Variable ( posvar, firstpoint, raf, lb, ub, varid, torestore ) |
---|
| 44 | !--------------------------------------------------------------------------------------------------- |
---|
| 45 | integer, dimension(:), intent(in) :: posvar !< position of the variable on the cell |
---|
| 46 | !! (1 for the border of the edge, 2 for the center) |
---|
| 47 | integer, dimension(:), intent(in) :: firstpoint !< index of the first point in the real domain |
---|
| 48 | character(1), dimension(:), intent(in) :: raf !< Array indicating the type of dimension (space or not) |
---|
| 49 | !! for each of them |
---|
| 50 | integer, dimension(:), intent(in) :: lb !< Lower bounds of the array |
---|
| 51 | integer, dimension(:), intent(in) :: ub !< Upper bounds of the array |
---|
| 52 | integer, intent(out) :: varid !< Id number of the newly created variable |
---|
| 53 | logical, optional, intent(in) :: torestore !< Indicates if the array restore is used |
---|
| 54 | !--------------------------------------------------------------------------------------------------- |
---|
| 55 | type(Agrif_Variables_List), pointer :: new_varlist |
---|
| 56 | type(Agrif_Variable), pointer :: var |
---|
| 57 | integer :: nbdim, i |
---|
| 58 | logical :: restore |
---|
| 59 | |
---|
| 60 | restore = .FALSE. |
---|
| 61 | if ( Agrif_Mygrid % ngridstep /= 0 ) then |
---|
| 62 | if (present(torestore)) restore = torestore |
---|
| 63 | endif |
---|
| 64 | ! |
---|
| 65 | nbdim = SIZE(posvar) |
---|
| 66 | ! |
---|
| 67 | allocate(new_varlist) |
---|
| 68 | allocate(new_varlist % var) |
---|
| 69 | |
---|
| 70 | var => new_varlist % var |
---|
| 71 | |
---|
| 72 | allocate(var % posvar(nbdim)) |
---|
| 73 | allocate(var % interptab(nbdim)) |
---|
| 74 | allocate(var % coords(nbdim)) |
---|
| 75 | ! |
---|
| 76 | var % nbdim = nbdim |
---|
| 77 | var % interptab = raf(1:nbdim) |
---|
| 78 | var % posvar = posvar(1:nbdim) |
---|
| 79 | var % point(1:nbdim) = firstpoint(1:nbdim) |
---|
| 80 | var % restore = restore |
---|
| 81 | ! |
---|
| 82 | do i = 1,nbdim |
---|
| 83 | select case( raf(i) ) |
---|
| 84 | case('x') ; var % coords(i) = 1 |
---|
| 85 | case('y') ; var % coords(i) = 2 |
---|
| 86 | case('z') ; var % coords(i) = 3 |
---|
| 87 | case('N') ; var % coords(i) = 0 |
---|
| 88 | case default ; var % coords(i) = 0 |
---|
| 89 | end select |
---|
| 90 | enddo |
---|
| 91 | ! |
---|
| 92 | var % lb(1:nbdim) = lb(1:nbdim) |
---|
| 93 | var % ub(1:nbdim) = ub(1:nbdim) |
---|
| 94 | |
---|
| 95 | |
---|
| 96 | if ( restore ) then |
---|
| 97 | select case(nbdim) |
---|
| 98 | case(1) |
---|
| 99 | allocate(var % Restore1D(lb(1):ub(1))) |
---|
| 100 | var % Restore1D = 0 |
---|
| 101 | case(2) |
---|
| 102 | allocate(var % Restore2D(lb(1):ub(1), & |
---|
| 103 | lb(2):ub(2))) |
---|
| 104 | var % Restore2D = 0 |
---|
| 105 | case(3) |
---|
| 106 | allocate(var % Restore3D(lb(1):ub(1), & |
---|
| 107 | lb(2):ub(2), & |
---|
| 108 | lb(3):ub(3))) |
---|
| 109 | var % Restore3D = 0 |
---|
| 110 | case(4) |
---|
| 111 | allocate(var % Restore4D(lb(1):ub(1), & |
---|
| 112 | lb(2):ub(2), & |
---|
| 113 | lb(3):ub(3), & |
---|
| 114 | lb(4):ub(4))) |
---|
| 115 | var % Restore4D = 0 |
---|
| 116 | case(5) |
---|
| 117 | allocate(var % Restore5D(lb(1):ub(1), & |
---|
| 118 | lb(2):ub(2), & |
---|
| 119 | lb(3):ub(3), & |
---|
| 120 | lb(4):ub(4), & |
---|
| 121 | lb(5):ub(5))) |
---|
| 122 | var % Restore5D = 0 |
---|
| 123 | end select |
---|
| 124 | endif |
---|
| 125 | |
---|
| 126 | new_varlist % next => Agrif_Curgrid % variables |
---|
| 127 | |
---|
| 128 | Agrif_Curgrid % variables => new_varlist |
---|
| 129 | Agrif_Curgrid % Nbvariables = Agrif_Curgrid % Nbvariables + 1 |
---|
| 130 | |
---|
| 131 | varid = Agrif_Curgrid % Nbvariables |
---|
| 132 | |
---|
| 133 | var % parent_var => Agrif_Search_Variable(Agrif_Curgrid % parent, Agrif_Curgrid % nbvariables) |
---|
| 134 | var % root_var => Agrif_Search_Variable(Agrif_Mygrid, Agrif_Curgrid % nbvariables) |
---|
| 135 | |
---|
| 136 | #if defined AGRIF_MPI |
---|
| 137 | call Agrif_get_var_global_bounds(var,var%lubglob(1:nbdim,:),nbdim) |
---|
| 138 | #else |
---|
| 139 | var%lubglob(1:nbdim,1) = var % lb(1:nbdim) |
---|
| 140 | var%lubglob(1:nbdim,2) = var % ub(1:nbdim) |
---|
| 141 | #endif |
---|
| 142 | !--------------------------------------------------------------------------------------------------- |
---|
| 143 | end subroutine Agrif_Declare_Variable |
---|
| 144 | !=================================================================================================== |
---|
| 145 | ! |
---|
| 146 | !=================================================================================================== |
---|
| 147 | ! function Agrif_Search_Variable |
---|
| 148 | ! |
---|
| 149 | !> \cond Returns a pointer to the variable varid for the grid grid. |
---|
| 150 | !--------------------------------------------------------------------------------------------------- |
---|
| 151 | function Agrif_Search_Variable ( grid, varid ) result(outvar) |
---|
| 152 | !--------------------------------------------------------------------------------------------------- |
---|
| 153 | type(Agrif_Grid), pointer :: grid !< Pointer on the current grid. |
---|
| 154 | integer, intent(in) :: varid !< ID number of the variable we are looking for. |
---|
| 155 | ! |
---|
| 156 | type(Agrif_Variable), pointer :: outvar |
---|
| 157 | type(Agrif_Variables_List), pointer :: parcours |
---|
| 158 | integer :: nb, varidinv |
---|
| 159 | ! |
---|
| 160 | if ( .not.associated(grid) ) then |
---|
| 161 | outvar => NULL() |
---|
| 162 | return |
---|
| 163 | endif |
---|
| 164 | ! |
---|
| 165 | parcours => grid % variables |
---|
| 166 | |
---|
| 167 | if (.not. associated(parcours)) then ! can occur on the grand mother grid |
---|
| 168 | outvar => NULL() ! during the first call by agrif_mygrid |
---|
| 169 | return |
---|
| 170 | endif |
---|
| 171 | |
---|
| 172 | varidinv = 1 + grid % nbvariables - varid |
---|
| 173 | |
---|
| 174 | do nb = 1,varidinv-1 |
---|
| 175 | parcours => parcours % next |
---|
| 176 | enddo |
---|
| 177 | |
---|
| 178 | outvar => parcours % var |
---|
| 179 | !--------------------------------------------------------------------------------------------------- |
---|
| 180 | end function Agrif_Search_variable |
---|
| 181 | !=================================================================================================== |
---|
| 182 | ! |
---|
| 183 | !!\endcond |
---|
| 184 | end module Agrif_User_Variables |
---|