1 | ! |
---|
2 | ! $Id$ |
---|
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_Boundary. |
---|
25 | !> |
---|
26 | !> Contains subroutines to calculate the boundary conditions on the child grids from their |
---|
27 | !> parent grids. |
---|
28 | ! |
---|
29 | module Agrif_Boundary |
---|
30 | ! |
---|
31 | use Agrif_Interpolation |
---|
32 | ! |
---|
33 | implicit none |
---|
34 | ! |
---|
35 | contains |
---|
36 | ! |
---|
37 | !=================================================================================================== |
---|
38 | ! subroutine Agrif_CorrectVariable |
---|
39 | ! |
---|
40 | !> subroutine to calculate the boundary conditions on a fine grid |
---|
41 | !--------------------------------------------------------------------------------------------------- |
---|
42 | subroutine Agrif_CorrectVariable ( parent, child, pweight, weight, procname ) |
---|
43 | !--------------------------------------------------------------------------------------------------- |
---|
44 | type(Agrif_Variable), pointer :: parent !< Variable on the parent grid |
---|
45 | type(Agrif_Variable), pointer :: child !< Variable on the child grid |
---|
46 | logical :: pweight !< Indicates if weight is used for the time interpolation |
---|
47 | real :: weight !< Coefficient for the time interpolation |
---|
48 | procedure() :: procname !< Data recovery procedure |
---|
49 | ! |
---|
50 | type(Agrif_Grid) , pointer :: Agrif_Child_Gr, Agrif_Parent_Gr |
---|
51 | type(Agrif_Variable), pointer :: root_var ! Variable on the root grid |
---|
52 | integer :: nbdim ! Number of dimensions of the grid variable |
---|
53 | integer :: n |
---|
54 | integer, dimension(6) :: lb_child ! Index of the first point inside the domain for |
---|
55 | ! the child grid variable |
---|
56 | integer, dimension(6) :: lb_parent ! Index of the first point inside the domain for |
---|
57 | ! the parent grid variable |
---|
58 | integer, dimension(6) :: ub_child ! Upper bound on the child grid |
---|
59 | integer, dimension(6) :: nb_child ! Number of cells for child |
---|
60 | integer, dimension(6) :: posvartab_child ! Position of the variable on the cell |
---|
61 | integer, dimension(6) :: loctab_child ! Indicates if the child grid has a common border |
---|
62 | ! with the root grid |
---|
63 | real, dimension(6) :: s_child, s_parent ! Positions of the parent and child grids |
---|
64 | real, dimension(6) :: ds_child, ds_parent ! Space steps of the parent and child grids |
---|
65 | ! |
---|
66 | call PreProcessToInterpOrUpdate( parent, child, & |
---|
67 | nb_child, ub_child, & |
---|
68 | lb_child, lb_parent, & |
---|
69 | s_child, s_parent, & |
---|
70 | ds_child, ds_parent, nbdim, interp=.true.) |
---|
71 | root_var => child % root_var |
---|
72 | Agrif_Child_Gr => Agrif_Curgrid |
---|
73 | Agrif_Parent_Gr => Agrif_Curgrid % parent |
---|
74 | ! |
---|
75 | loctab_child(1:nbdim) = 0 |
---|
76 | posvartab_child(1:nbdim) = root_var % posvar(1:nbdim) |
---|
77 | ! |
---|
78 | do n = 1,nbdim |
---|
79 | ! |
---|
80 | select case(root_var % interptab(n)) |
---|
81 | ! |
---|
82 | case('x') ! x DIMENSION |
---|
83 | ! |
---|
84 | if (Agrif_Curgrid % NearRootBorder(1)) loctab_child(n) = -1 |
---|
85 | if (Agrif_Curgrid % DistantRootBorder(1)) loctab_child(n) = -2 |
---|
86 | if ((Agrif_Curgrid % NearRootBorder(1)) .AND. & |
---|
87 | (Agrif_Curgrid % DistantRootBorder(1))) loctab_child(n) = -3 |
---|
88 | ! |
---|
89 | case('y') ! y DIMENSION |
---|
90 | ! |
---|
91 | if (Agrif_Curgrid % NearRootBorder(2)) loctab_child(n) = -1 |
---|
92 | if (Agrif_Curgrid % DistantRootBorder(2)) loctab_child(n) = -2 |
---|
93 | if ((Agrif_Curgrid % NearRootBorder(2)) .AND. & |
---|
94 | (Agrif_Curgrid % DistantRootBorder(2))) loctab_child(n) = -3 |
---|
95 | ! |
---|
96 | case('z') ! z DIMENSION |
---|
97 | ! |
---|
98 | if (Agrif_Curgrid % NearRootBorder(3)) loctab_child(n) = -1 |
---|
99 | if (Agrif_Curgrid % DistantRootBorder(3)) loctab_child(n) = -2 |
---|
100 | if ((Agrif_Curgrid % NearRootBorder(3)) .AND. & |
---|
101 | (Agrif_Curgrid % DistantRootBorder(3))) loctab_child(n) = -3 |
---|
102 | ! |
---|
103 | case('N') ! No space DIMENSION |
---|
104 | ! |
---|
105 | posvartab_child(n) = 1 |
---|
106 | loctab_child(n) = -3 |
---|
107 | ! |
---|
108 | end select |
---|
109 | ! |
---|
110 | enddo |
---|
111 | ! |
---|
112 | call Agrif_Correctnd(parent, child, pweight, weight, & |
---|
113 | lb_child(1:nbdim), lb_parent(1:nbdim), & |
---|
114 | nb_child(1:nbdim), posvartab_child(1:nbdim), & |
---|
115 | loctab_child(1:nbdim), & |
---|
116 | s_child(1:nbdim), s_parent(1:nbdim), & |
---|
117 | ds_child(1:nbdim),ds_parent(1:nbdim), nbdim, procname ) |
---|
118 | !--------------------------------------------------------------------------------------------------- |
---|
119 | end subroutine Agrif_CorrectVariable |
---|
120 | !=================================================================================================== |
---|
121 | ! |
---|
122 | !=================================================================================================== |
---|
123 | ! subroutine Agrif_Correctnd |
---|
124 | ! |
---|
125 | !> calculates the boundary conditions for a nD grid variable on a fine grid by using |
---|
126 | !> a space and time interpolations; it is called by the #Agrif_CorrectVariable procedure |
---|
127 | !--------------------------------------------------------------------------------------------------- |
---|
128 | subroutine Agrif_Correctnd ( parent, child, pweight, weight, & |
---|
129 | pttab_child, pttab_Parent, & |
---|
130 | nbtab_Child, posvartab_Child, loctab_Child, & |
---|
131 | s_Child, s_Parent, ds_Child, ds_Parent, & |
---|
132 | nbdim, procname ) |
---|
133 | !--------------------------------------------------------------------------------------------------- |
---|
134 | #if defined AGRIF_MPI |
---|
135 | include 'mpif.h' |
---|
136 | #endif |
---|
137 | ! |
---|
138 | TYPE(Agrif_Variable), pointer :: parent !< Variable on the parent grid |
---|
139 | TYPE(Agrif_Variable), pointer :: child !< Variable on the child grid |
---|
140 | LOGICAL :: pweight !< Indicates if weight is used for the temporal interpolation |
---|
141 | REAL :: weight !< Coefficient for the temporal interpolation |
---|
142 | INTEGER, DIMENSION(nbdim) :: pttab_child !< Index of the first point inside the domain for the parent grid variable |
---|
143 | INTEGER, DIMENSION(nbdim) :: pttab_Parent !< Index of the first point inside the domain for the child grid variable |
---|
144 | INTEGER, DIMENSION(nbdim) :: nbtab_Child !< Number of cells of the child grid |
---|
145 | INTEGER, DIMENSION(nbdim) :: posvartab_Child !< Position of the grid variable (1 or 2) |
---|
146 | INTEGER, DIMENSION(nbdim) :: loctab_Child !< Indicates if the child grid has a common border with the root grid |
---|
147 | REAL , DIMENSION(nbdim) :: s_Child, s_Parent !< Positions of the parent and child grids |
---|
148 | REAL , DIMENSION(nbdim) :: ds_Child, ds_Parent !< Space steps of the parent and child grids |
---|
149 | INTEGER :: nbdim !< Number of dimensions of the grid variable |
---|
150 | procedure() :: procname !< Data recovery procedure |
---|
151 | ! |
---|
152 | INTEGER,DIMENSION(6) :: type_interp ! Type of interpolation (linear, spline,...) |
---|
153 | INTEGER,DIMENSION(6,6) :: type_interp_bc ! Type of interpolation (linear, spline,...) |
---|
154 | INTEGER,DIMENSION(nbdim,2,2) :: childarray |
---|
155 | INTEGER,DIMENSION(nbdim,2) :: lubglob |
---|
156 | INTEGER :: kindex ! Index used for safeguard and time interpolation |
---|
157 | INTEGER,DIMENSION(nbdim,2,2) :: indtab ! Arrays indicating the limits of the child |
---|
158 | INTEGER,DIMENSION(nbdim,2,2) :: indtruetab ! grid variable where boundary conditions are |
---|
159 | INTEGER,DIMENSION(nbdim,2,2,nbdim) :: ptres,ptres2 ! calculated |
---|
160 | INTEGER,DIMENSION(nbdim) :: coords |
---|
161 | INTEGER :: i, nb, ndir |
---|
162 | INTEGER :: n, sizetab |
---|
163 | INTEGER :: ibeg, iend |
---|
164 | INTEGER :: i1,i2,j1,j2,k1,k2,l1,l2,m1,m2,n1,n2 |
---|
165 | REAL :: c1t,c2t ! Coefficients for the time interpolation (c2t=1-c1t) |
---|
166 | #if defined AGRIF_MPI |
---|
167 | ! |
---|
168 | INTEGER, DIMENSION(nbdim) :: lower, upper |
---|
169 | INTEGER, DIMENSION(nbdim) :: ltab, utab |
---|
170 | ! |
---|
171 | #endif |
---|
172 | ! |
---|
173 | type_interp_bc = child % root_var % type_interp_bc |
---|
174 | coords = child % root_var % coords |
---|
175 | ! |
---|
176 | ibeg = child % bcinf |
---|
177 | iend = child % bcsup |
---|
178 | ! |
---|
179 | indtab(1:nbdim,2,1) = pttab_child(1:nbdim) + nbtab_child(1:nbdim) + ibeg |
---|
180 | indtab(1:nbdim,2,2) = indtab(1:nbdim,2,1) + ( iend - ibeg ) |
---|
181 | |
---|
182 | indtab(1:nbdim,1,1) = pttab_child(1:nbdim) - iend |
---|
183 | indtab(1:nbdim,1,2) = pttab_child(1:nbdim) - ibeg |
---|
184 | |
---|
185 | WHERE (posvartab_child(1:nbdim) == 2) |
---|
186 | indtab(1:nbdim,1,1) = indtab(1:nbdim,1,1) - 1 |
---|
187 | indtab(1:nbdim,1,2) = indtab(1:nbdim,1,2) - 1 |
---|
188 | END WHERE |
---|
189 | ! |
---|
190 | call Agrif_get_var_global_bounds(child,lubglob,nbdim,parent) |
---|
191 | ! |
---|
192 | indtruetab(1:nbdim,1,1) = max(indtab(1:nbdim,1,1), lubglob(1:nbdim,1)) |
---|
193 | indtruetab(1:nbdim,1,2) = max(indtab(1:nbdim,1,2), lubglob(1:nbdim,1)) |
---|
194 | indtruetab(1:nbdim,2,1) = min(indtab(1:nbdim,2,1), lubglob(1:nbdim,2)) |
---|
195 | indtruetab(1:nbdim,2,2) = min(indtab(1:nbdim,2,2), lubglob(1:nbdim,2)) |
---|
196 | ! |
---|
197 | do nb = 1,nbdim |
---|
198 | do ndir = 1,2 |
---|
199 | ! |
---|
200 | if (loctab_child(nb) /= (-ndir) .AND. loctab_child(nb) /= -3) then |
---|
201 | ! |
---|
202 | do n = 1,2 |
---|
203 | ptres(nb,n,ndir,nb) = indtruetab(nb,ndir,n) |
---|
204 | enddo |
---|
205 | ! |
---|
206 | do i = 1,nbdim |
---|
207 | ! |
---|
208 | if (i /= nb) then |
---|
209 | ! |
---|
210 | if (loctab_child(i) == -1 .OR. loctab_child(i) == -3) then |
---|
211 | ptres(i,1,ndir,nb) = pttab_child(i) |
---|
212 | else |
---|
213 | ptres(i,1,ndir,nb) = indtruetab(i,1,1) |
---|
214 | endif |
---|
215 | if (loctab_child(i) == -2 .OR. loctab_child(i) == -3) then |
---|
216 | if (posvartab_child(i) == 1) then |
---|
217 | ptres(i,2,ndir,nb) = pttab_child(i) + nbtab_child(i) |
---|
218 | else |
---|
219 | ptres(i,2,ndir,nb) = pttab_child(i) + nbtab_child(i) - 1 |
---|
220 | endif |
---|
221 | else |
---|
222 | ptres(i,2,ndir,nb) = indtruetab(i,2,2) |
---|
223 | endif |
---|
224 | ! |
---|
225 | endif |
---|
226 | ! |
---|
227 | enddo |
---|
228 | |
---|
229 | ! |
---|
230 | #if defined AGRIF_MPI |
---|
231 | call Agrif_get_var_bounds_array(child,lower,upper,nbdim,parent) |
---|
232 | |
---|
233 | do i = 1,nbdim |
---|
234 | ! |
---|
235 | Call Agrif_GetLocalBoundaries(ptres(i,1,ndir,nb), ptres(i,2,ndir,nb), & |
---|
236 | coords(i), lower(i), upper(i), ltab(i), utab(i) ) |
---|
237 | ptres2(i,1,ndir,nb) = max(ltab(i),lower(i)) |
---|
238 | ptres2(i,2,ndir,nb) = min(utab(i),upper(i)) |
---|
239 | if ((i == nb) .AND. (ndir == 1)) then |
---|
240 | ptres2(i,2,ndir,nb) = max(utab(i),lower(i)) |
---|
241 | elseif ((i == nb) .AND. (ndir == 2)) then |
---|
242 | ptres2(i,1,ndir,nb) = min(ltab(i),upper(i)) |
---|
243 | endif |
---|
244 | ! |
---|
245 | enddo |
---|
246 | #else |
---|
247 | ptres2(:,:,ndir,nb) = ptres(:,:,ndir,nb) |
---|
248 | #endif |
---|
249 | endif |
---|
250 | ! |
---|
251 | enddo ! ndir = 1,2 |
---|
252 | enddo ! nb = 1,nbdim |
---|
253 | ! |
---|
254 | if ( child % interpIndex /= Agrif_Curgrid % parent % ngridstep .OR. & |
---|
255 | child % Interpolationshouldbemade ) then |
---|
256 | ! |
---|
257 | ! Space interpolation |
---|
258 | ! |
---|
259 | kindex = 1 |
---|
260 | ! |
---|
261 | do nb = 1,nbdim |
---|
262 | |
---|
263 | type_interp = type_interp_bc(nb,:) |
---|
264 | |
---|
265 | do ndir = 1,2 |
---|
266 | ! |
---|
267 | if (loctab_child(nb) /= (-ndir) .AND. loctab_child(nb) /= -3) then |
---|
268 | ! |
---|
269 | call Agrif_InterpnD(type_interp, parent, child, & |
---|
270 | ptres(1:nbdim,1,ndir,nb), & |
---|
271 | ptres(1:nbdim,2,ndir,nb), & |
---|
272 | pttab_child(1:nbdim), & |
---|
273 | pttab_Parent(1:nbdim), & |
---|
274 | s_Child(1:nbdim), s_Parent(1:nbdim), & |
---|
275 | ds_Child(1:nbdim),ds_Parent(1:nbdim), & |
---|
276 | NULL(), .FALSE., nbdim, & |
---|
277 | childarray, & |
---|
278 | child%memberin(nb,ndir), .TRUE., procname, coords(nb), ndir) |
---|
279 | |
---|
280 | child % childarray(1:nbdim,:,:,nb,ndir) = childarray |
---|
281 | |
---|
282 | if (.not. child%interpolationshouldbemade) then |
---|
283 | ! |
---|
284 | ! Safeguard of the values of the grid variable (at times n and n+1 on the parent grid) |
---|
285 | ! |
---|
286 | sizetab = 1 |
---|
287 | do i = 1,nbdim |
---|
288 | sizetab = sizetab * (ptres2(i,2,ndir,nb)-ptres2(i,1,ndir,nb)+1) |
---|
289 | enddo |
---|
290 | |
---|
291 | call saveAfterInterp(child,ptres2(:,:,ndir,nb),kindex,sizetab,nbdim) |
---|
292 | ! |
---|
293 | endif |
---|
294 | ! |
---|
295 | endif |
---|
296 | ! |
---|
297 | enddo ! ndir = 1,2 |
---|
298 | enddo ! nb = 1,nbdim |
---|
299 | ! |
---|
300 | child % interpIndex = Agrif_Curgrid % parent % ngridstep |
---|
301 | ! |
---|
302 | endif |
---|
303 | ! |
---|
304 | if (.not. child%interpolationshouldbemade) then |
---|
305 | ! |
---|
306 | ! Calculation of the coefficients c1t and c2t for the temporary interpolation |
---|
307 | ! |
---|
308 | if (pweight) then |
---|
309 | c1t = weight |
---|
310 | else |
---|
311 | c1t = (REAL(Agrif_Nbstepint()) + 1.) / Agrif_Rhot() |
---|
312 | endif |
---|
313 | c2t = 1. - c1t |
---|
314 | ! |
---|
315 | ! Time interpolation |
---|
316 | ! |
---|
317 | kindex = 1 |
---|
318 | ! |
---|
319 | do nb = 1,nbdim |
---|
320 | do ndir = 1,2 |
---|
321 | if (loctab_child(nb) /= (-ndir) .AND. loctab_child(nb) /= -3) then |
---|
322 | Call timeInterpolation(child,ptres2(:,:,ndir,nb),kindex,c1t,c2t,nbdim) |
---|
323 | endif |
---|
324 | enddo |
---|
325 | enddo |
---|
326 | ! |
---|
327 | endif |
---|
328 | ! |
---|
329 | do nb = 1,nbdim |
---|
330 | do ndir = 1,2 |
---|
331 | if ( (loctab_child(nb) /= (-ndir)) .AND. (loctab_child(nb) /= -3) .AND. child%memberin(nb,ndir) ) then |
---|
332 | select case(nbdim) |
---|
333 | case(1) |
---|
334 | i1 = child % childarray(1,1,2,nb,ndir) |
---|
335 | i2 = child % childarray(1,2,2,nb,ndir) |
---|
336 | |
---|
337 | call procname(parray1(i1:i2), & |
---|
338 | i1,i2, .FALSE.,coords(nb),ndir) |
---|
339 | case(2) |
---|
340 | i1 = child % childarray(1,1,2,nb,ndir) |
---|
341 | i2 = child % childarray(1,2,2,nb,ndir) |
---|
342 | j1 = child % childarray(2,1,2,nb,ndir) |
---|
343 | j2 = child % childarray(2,2,2,nb,ndir) |
---|
344 | |
---|
345 | call procname(parray2(i1:i2,j1:j2), & |
---|
346 | i1,i2,j1,j2, .FALSE.,coords(nb),ndir) |
---|
347 | case(3) |
---|
348 | i1 = child % childarray(1,1,2,nb,ndir) |
---|
349 | i2 = child % childarray(1,2,2,nb,ndir) |
---|
350 | j1 = child % childarray(2,1,2,nb,ndir) |
---|
351 | j2 = child % childarray(2,2,2,nb,ndir) |
---|
352 | k1 = child % childarray(3,1,2,nb,ndir) |
---|
353 | k2 = child % childarray(3,2,2,nb,ndir) |
---|
354 | |
---|
355 | call procname(parray3(i1:i2,j1:j2,k1:k2), & |
---|
356 | i1,i2,j1,j2,k1,k2, .FALSE.,coords(nb),ndir) |
---|
357 | case(4) |
---|
358 | i1 = child % childarray(1,1,2,nb,ndir) |
---|
359 | i2 = child % childarray(1,2,2,nb,ndir) |
---|
360 | j1 = child % childarray(2,1,2,nb,ndir) |
---|
361 | j2 = child % childarray(2,2,2,nb,ndir) |
---|
362 | k1 = child % childarray(3,1,2,nb,ndir) |
---|
363 | k2 = child % childarray(3,2,2,nb,ndir) |
---|
364 | l1 = child % childarray(4,1,2,nb,ndir) |
---|
365 | l2 = child % childarray(4,2,2,nb,ndir) |
---|
366 | |
---|
367 | call procname(parray4(i1:i2,j1:j2,k1:k2,l1:l2), & |
---|
368 | i1,i2,j1,j2,k1,k2,l1,l2, .FALSE.,coords(nb),ndir) |
---|
369 | case(5) |
---|
370 | i1 = child % childarray(1,1,2,nb,ndir) |
---|
371 | i2 = child % childarray(1,2,2,nb,ndir) |
---|
372 | j1 = child % childarray(2,1,2,nb,ndir) |
---|
373 | j2 = child % childarray(2,2,2,nb,ndir) |
---|
374 | k1 = child % childarray(3,1,2,nb,ndir) |
---|
375 | k2 = child % childarray(3,2,2,nb,ndir) |
---|
376 | l1 = child % childarray(4,1,2,nb,ndir) |
---|
377 | l2 = child % childarray(4,2,2,nb,ndir) |
---|
378 | m1 = child % childarray(5,1,2,nb,ndir) |
---|
379 | m2 = child % childarray(5,2,2,nb,ndir) |
---|
380 | |
---|
381 | call procname(parray5(i1:i2,j1:j2,k1:k2,l1:l2,m1:m2), & |
---|
382 | i1,i2,j1,j2,k1,k2,l1,l2,m1,m2, .FALSE.,coords(nb),ndir) |
---|
383 | case(6) |
---|
384 | i1 = child % childarray(1,1,2,nb,ndir) |
---|
385 | i2 = child % childarray(1,2,2,nb,ndir) |
---|
386 | j1 = child % childarray(2,1,2,nb,ndir) |
---|
387 | j2 = child % childarray(2,2,2,nb,ndir) |
---|
388 | k1 = child % childarray(3,1,2,nb,ndir) |
---|
389 | k2 = child % childarray(3,2,2,nb,ndir) |
---|
390 | l1 = child % childarray(4,1,2,nb,ndir) |
---|
391 | l2 = child % childarray(4,2,2,nb,ndir) |
---|
392 | m1 = child % childarray(5,1,2,nb,ndir) |
---|
393 | m2 = child % childarray(5,2,2,nb,ndir) |
---|
394 | n1 = child % childarray(6,1,2,nb,ndir) |
---|
395 | n2 = child % childarray(6,2,2,nb,ndir) |
---|
396 | |
---|
397 | call procname(parray6(i1:i2,j1:j2,k1:k2,l1:l2,m1:m2,n1:n2), & |
---|
398 | i1,i2,j1,j2,k1,k2,l1,l2,m1,m2,n1,n2, .FALSE.,coords(nb),ndir) |
---|
399 | end select |
---|
400 | endif |
---|
401 | enddo |
---|
402 | enddo |
---|
403 | !--------------------------------------------------------------------------------------------------- |
---|
404 | end subroutine Agrif_Correctnd |
---|
405 | !=================================================================================================== |
---|
406 | ! |
---|
407 | !=================================================================================================== |
---|
408 | ! subroutine saveAfterInterp |
---|
409 | ! |
---|
410 | !> saves the values of the grid variable on the fine grid after the space interpolation |
---|
411 | !--------------------------------------------------------------------------------------------------- |
---|
412 | subroutine saveAfterInterp ( child_var, bounds, kindex, newsize, nbdim ) |
---|
413 | !--------------------------------------------------------------------------------------------------- |
---|
414 | TYPE (Agrif_Variable), INTENT(inout) :: child_var !< The fine grid variable |
---|
415 | INTEGER, DIMENSION(nbdim,2), INTENT(in) :: bounds |
---|
416 | INTEGER, INTENT(inout) :: kindex !< Index indicating where this safeguard |
---|
417 | !< is done on the fine grid |
---|
418 | INTEGER, INTENT(in) :: newsize |
---|
419 | INTEGER, INTENT(in) :: nbdim |
---|
420 | ! |
---|
421 | INTEGER :: ir,jr,kr,lr,mr,nr |
---|
422 | ! |
---|
423 | ! Allocation of the array oldvalues2d |
---|
424 | ! |
---|
425 | if (newsize .LE. 0) return |
---|
426 | ! |
---|
427 | Call Agrif_Checksize(child_var,kindex+newsize) |
---|
428 | |
---|
429 | if (child_var % interpIndex /= Agrif_Curgrid % parent % ngridstep ) then |
---|
430 | child_var % oldvalues2d(1,kindex:kindex+newsize-1) = & |
---|
431 | child_var % oldvalues2d(2,kindex:kindex+newsize-1) |
---|
432 | endif |
---|
433 | |
---|
434 | SELECT CASE (nbdim) |
---|
435 | CASE (1) |
---|
436 | !CDIR ALTCODE |
---|
437 | do ir = bounds(1,1), bounds(1,2) |
---|
438 | child_var % oldvalues2d(2,kindex) = parray1(ir) |
---|
439 | kindex = kindex + 1 |
---|
440 | enddo |
---|
441 | ! |
---|
442 | CASE (2) |
---|
443 | do jr = bounds(2,1),bounds(2,2) |
---|
444 | !CDIR ALTCODE |
---|
445 | do ir = bounds(1,1),bounds(1,2) |
---|
446 | child_var % oldvalues2d(2,kindex) = parray2(ir,jr) |
---|
447 | kindex = kindex + 1 |
---|
448 | enddo |
---|
449 | enddo |
---|
450 | ! |
---|
451 | CASE (3) |
---|
452 | do kr = bounds(3,1),bounds(3,2) |
---|
453 | do jr = bounds(2,1),bounds(2,2) |
---|
454 | !CDIR ALTCODE |
---|
455 | do ir = bounds(1,1),bounds(1,2) |
---|
456 | child_var % oldvalues2d(2,kindex) = parray3(ir,jr,kr) |
---|
457 | kindex = kindex + 1 |
---|
458 | enddo |
---|
459 | enddo |
---|
460 | enddo |
---|
461 | ! |
---|
462 | CASE (4) |
---|
463 | do lr = bounds(4,1),bounds(4,2) |
---|
464 | do kr = bounds(3,1),bounds(3,2) |
---|
465 | do jr = bounds(2,1),bounds(2,2) |
---|
466 | !CDIR ALTCODE |
---|
467 | do ir = bounds(1,1),bounds(1,2) |
---|
468 | child_var % oldvalues2d(2,kindex) = parray4(ir,jr,kr,lr) |
---|
469 | kindex = kindex + 1 |
---|
470 | enddo |
---|
471 | enddo |
---|
472 | enddo |
---|
473 | enddo |
---|
474 | ! |
---|
475 | CASE (5) |
---|
476 | do mr = bounds(5,1),bounds(5,2) |
---|
477 | do lr = bounds(4,1),bounds(4,2) |
---|
478 | do kr = bounds(3,1),bounds(3,2) |
---|
479 | do jr = bounds(2,1),bounds(2,2) |
---|
480 | !CDIR ALTCODE |
---|
481 | do ir = bounds(1,1),bounds(1,2) |
---|
482 | child_var % oldvalues2d(2,kindex) = parray5(ir,jr,kr,lr,mr) |
---|
483 | kindex = kindex + 1 |
---|
484 | enddo |
---|
485 | enddo |
---|
486 | enddo |
---|
487 | enddo |
---|
488 | enddo |
---|
489 | ! |
---|
490 | CASE (6) |
---|
491 | do nr = bounds(6,1),bounds(6,2) |
---|
492 | do mr = bounds(5,1),bounds(5,2) |
---|
493 | do lr = bounds(4,1),bounds(4,2) |
---|
494 | do kr = bounds(3,1),bounds(3,2) |
---|
495 | do jr = bounds(2,1),bounds(2,2) |
---|
496 | !CDIR ALTCODE |
---|
497 | do ir = bounds(1,1),bounds(1,2) |
---|
498 | child_var % oldvalues2d(2,kindex) = parray6(ir,jr,kr,lr,mr,nr) |
---|
499 | kindex = kindex + 1 |
---|
500 | enddo |
---|
501 | enddo |
---|
502 | enddo |
---|
503 | enddo |
---|
504 | enddo |
---|
505 | enddo |
---|
506 | END SELECT |
---|
507 | !--------------------------------------------------------------------------------------------------- |
---|
508 | end subroutine saveAfterInterp |
---|
509 | !=================================================================================================== |
---|
510 | ! |
---|
511 | !=================================================================================================== |
---|
512 | ! subroutine timeInterpolation |
---|
513 | ! |
---|
514 | !> subroutine for a linear time interpolation on the child grid |
---|
515 | !--------------------------------------------------------------------------------------------------- |
---|
516 | subroutine timeInterpolation ( child_var, bounds, kindex, c1t, c2t, nbdim ) |
---|
517 | !--------------------------------------------------------------------------------------------------- |
---|
518 | TYPE (Agrif_Variable) :: child_var !< The fine grid variable |
---|
519 | INTEGER, DIMENSION(nbdim,2) :: bounds |
---|
520 | INTEGER :: kindex !< Index indicating the values of the fine grid got |
---|
521 | !< before and after the space interpolation and |
---|
522 | !< used for the time interpolation |
---|
523 | REAL :: c1t, c2t !< Coefficients for the time interpolation (c2t=1-c1t) |
---|
524 | INTEGER :: nbdim |
---|
525 | ! |
---|
526 | INTEGER :: ir,jr,kr,lr,mr,nr |
---|
527 | ! |
---|
528 | SELECT CASE (nbdim) |
---|
529 | CASE (1) |
---|
530 | !CDIR ALTCODE |
---|
531 | do ir = bounds(1,1),bounds(1,2) |
---|
532 | parray1(ir) = c2t*child_var % oldvalues2d(1,kindex) + & |
---|
533 | c1t*child_var % oldvalues2d(2,kindex) |
---|
534 | kindex = kindex + 1 |
---|
535 | enddo |
---|
536 | ! |
---|
537 | CASE (2) |
---|
538 | do jr = bounds(2,1),bounds(2,2) |
---|
539 | !CDIR ALTCODE |
---|
540 | do ir = bounds(1,1),bounds(1,2) |
---|
541 | parray2(ir,jr) = c2t*child_var % oldvalues2d(1,kindex) + & |
---|
542 | c1t*child_var % oldvalues2d(2,kindex) |
---|
543 | kindex = kindex + 1 |
---|
544 | enddo |
---|
545 | enddo |
---|
546 | ! |
---|
547 | CASE (3) |
---|
548 | do kr = bounds(3,1),bounds(3,2) |
---|
549 | do jr = bounds(2,1),bounds(2,2) |
---|
550 | !CDIR ALTCODE |
---|
551 | do ir = bounds(1,1),bounds(1,2) |
---|
552 | parray3(ir,jr,kr) = c2t*child_var % oldvalues2d(1,kindex) + & |
---|
553 | c1t*child_var % oldvalues2d(2,kindex) |
---|
554 | kindex = kindex + 1 |
---|
555 | enddo |
---|
556 | enddo |
---|
557 | enddo |
---|
558 | ! |
---|
559 | CASE (4) |
---|
560 | do lr = bounds(4,1),bounds(4,2) |
---|
561 | do kr = bounds(3,1),bounds(3,2) |
---|
562 | do jr = bounds(2,1),bounds(2,2) |
---|
563 | !CDIR ALTCODE |
---|
564 | do ir = bounds(1,1),bounds(1,2) |
---|
565 | parray4(ir,jr,kr,lr) = c2t*child_var % oldvalues2d(1,kindex) + & |
---|
566 | c1t*child_var % oldvalues2d(2,kindex) |
---|
567 | kindex = kindex + 1 |
---|
568 | enddo |
---|
569 | enddo |
---|
570 | enddo |
---|
571 | enddo |
---|
572 | ! |
---|
573 | CASE (5) |
---|
574 | do mr=bounds(5,1),bounds(5,2) |
---|
575 | do lr=bounds(4,1),bounds(4,2) |
---|
576 | do kr=bounds(3,1),bounds(3,2) |
---|
577 | do jr=bounds(2,1),bounds(2,2) |
---|
578 | !CDIR ALTCODE |
---|
579 | do ir=bounds(1,1),bounds(1,2) |
---|
580 | parray5(ir,jr,kr,lr,mr) = c2t*child_var % oldvalues2d(1,kindex) + & |
---|
581 | c1t*child_var % oldvalues2d(2,kindex) |
---|
582 | kindex = kindex + 1 |
---|
583 | enddo |
---|
584 | enddo |
---|
585 | enddo |
---|
586 | enddo |
---|
587 | enddo |
---|
588 | ! |
---|
589 | CASE (6) |
---|
590 | do nr=bounds(6,1),bounds(6,2) |
---|
591 | do mr=bounds(5,1),bounds(5,2) |
---|
592 | do lr=bounds(4,1),bounds(4,2) |
---|
593 | do kr=bounds(3,1),bounds(3,2) |
---|
594 | do jr=bounds(2,1),bounds(2,2) |
---|
595 | !CDIR ALTCODE |
---|
596 | do ir=bounds(1,1),bounds(1,2) |
---|
597 | parray6(ir,jr,kr,lr,mr,nr) = c2t*child_var % oldvalues2d(1,kindex) + & |
---|
598 | c1t*child_var % oldvalues2d(2,kindex) |
---|
599 | kindex = kindex + 1 |
---|
600 | enddo |
---|
601 | enddo |
---|
602 | enddo |
---|
603 | enddo |
---|
604 | enddo |
---|
605 | enddo |
---|
606 | END SELECT |
---|
607 | !--------------------------------------------------------------------------------------------------- |
---|
608 | end subroutine timeInterpolation |
---|
609 | !=================================================================================================== |
---|
610 | ! |
---|
611 | !=================================================================================================== |
---|
612 | ! subroutine Agrif_Checksize |
---|
613 | ! |
---|
614 | !> subroutine used in the saveAfterInterp procedure to allocate the oldvalues2d array |
---|
615 | !--------------------------------------------------------------------------------------------------- |
---|
616 | subroutine Agrif_Checksize ( child_var, newsize ) |
---|
617 | !--------------------------------------------------------------------------------------------------- |
---|
618 | TYPE (Agrif_Variable), INTENT(inout) :: child_var !< The fine grid variable |
---|
619 | INTEGER , INTENT(in) :: newsize !< Size of the domains where the boundary |
---|
620 | !< conditions are calculated |
---|
621 | ! |
---|
622 | REAL, DIMENSION(:,:), Allocatable :: tempoldvalues ! Temporary array |
---|
623 | ! |
---|
624 | if (.NOT. associated(child_var % oldvalues2d)) then |
---|
625 | ! |
---|
626 | allocate(child_var % oldvalues2d(2,newsize)) |
---|
627 | child_var % oldvalues2d = 0. |
---|
628 | ! |
---|
629 | else |
---|
630 | ! |
---|
631 | if (SIZE(child_var % oldvalues2d,2) < newsize) then |
---|
632 | ! |
---|
633 | allocate(tempoldvalues(2,SIZE(child_var % oldvalues2d,2))) |
---|
634 | tempoldvalues = child_var % oldvalues2d |
---|
635 | deallocate(child_var % oldvalues2d) |
---|
636 | allocate( child_var % oldvalues2d(2,newsize)) |
---|
637 | child_var % oldvalues2d = 0. |
---|
638 | child_var % oldvalues2d(:,1:SIZE(tempoldvalues,2)) = tempoldvalues(:,:) |
---|
639 | deallocate(tempoldvalues) |
---|
640 | ! |
---|
641 | endif |
---|
642 | ! |
---|
643 | endif |
---|
644 | !--------------------------------------------------------------------------------------------------- |
---|
645 | end subroutine Agrif_Checksize |
---|
646 | !=================================================================================================== |
---|
647 | ! |
---|
648 | end module Agrif_Boundary |
---|
649 | |
---|