1 | ! |
---|
2 | ! $Id$ |
---|
3 | ! |
---|
4 | C AGRIF (Adaptive Grid Refinement In Fortran) |
---|
5 | C |
---|
6 | C Copyright (C) 2003 Laurent Debreu (Laurent.Debreu@imag.fr) |
---|
7 | C Christophe Vouland (Christophe.Vouland@imag.fr) |
---|
8 | C |
---|
9 | C This program is free software; you can redistribute it and/or modify |
---|
10 | C it under the terms of the GNU General Public License as published by |
---|
11 | C the Free Software Foundation; either version 2 of the License, or |
---|
12 | C (at your option) any later version. |
---|
13 | C |
---|
14 | C This program is distributed in the hope that it will be useful, |
---|
15 | C but WITHOUT ANY WARRANTY; without even the implied warranty of |
---|
16 | C MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the |
---|
17 | C GNU General Public License for more details. |
---|
18 | C |
---|
19 | C You should have received a copy of the GNU General Public License |
---|
20 | C along with this program; if not, write to the Free Software |
---|
21 | C Foundation, Inc., 59 Temple Place- Suite 330, Boston, MA 02111-1307, USA. |
---|
22 | C |
---|
23 | C |
---|
24 | C |
---|
25 | CCC Module Agrif_Boundary |
---|
26 | C |
---|
27 | Module Agrif_Boundary |
---|
28 | C |
---|
29 | CCC Description: |
---|
30 | CCC Module to calculate the boundary conditions on the child grids from their |
---|
31 | CCC parent grids. |
---|
32 | C |
---|
33 | C Modules used: |
---|
34 | C |
---|
35 | Use Agrif_Interpolation |
---|
36 | C |
---|
37 | IMPLICIT NONE |
---|
38 | C |
---|
39 | CONTAINS |
---|
40 | C Define procedures contained in this module |
---|
41 | C |
---|
42 | C |
---|
43 | |
---|
44 | C |
---|
45 | C ************************************************************************** |
---|
46 | CCC Subroutine Agrif_Interp_bc_1d |
---|
47 | C ************************************************************************** |
---|
48 | C |
---|
49 | Subroutine Agrif_Interp_bc_1d(TypeInterp,parent,child,tab,deb,fin, |
---|
50 | & weight,pweight,procname) |
---|
51 | C |
---|
52 | CCC Description: |
---|
53 | CCC Subroutine to calculate the boundary conditions on a fine grid for a 1D |
---|
54 | CCC grid variable. |
---|
55 | C |
---|
56 | C Declarations: |
---|
57 | C |
---|
58 | |
---|
59 | C |
---|
60 | C Arguments |
---|
61 | External :: procname |
---|
62 | Optional :: procname |
---|
63 | INTEGER,DIMENSION(6,6) :: TypeInterp ! TYPE of interpolation (linear, |
---|
64 | ! lagrange, spline, ... ) |
---|
65 | TYPE(AGRIF_PVariable) :: parent ! Variable on the parent grid |
---|
66 | TYPE(AGRIF_PVariable) :: child ! Variable on the child grid |
---|
67 | TYPE(AGRIF_PVariable) :: childtemp ! Temporary variable on the child |
---|
68 | ! grid |
---|
69 | INTEGER :: deb,fin ! Positions where interpolations are |
---|
70 | ! done on the fine grid |
---|
71 | REAL, DIMENSION( |
---|
72 | & child%var%lb(1):child%var%ub(1) |
---|
73 | & ), Target :: tab ! Values of the grid variable |
---|
74 | LOGICAL :: pweight ! Indicates if weight is used for |
---|
75 | ! the temporal interpolation |
---|
76 | REAL :: weight ! Coefficient for the time |
---|
77 | ! interpolation |
---|
78 | C |
---|
79 | C |
---|
80 | C Definition of a temporary AGRIF_PVariable data TYPE representing the grid |
---|
81 | C variable. |
---|
82 | C |
---|
83 | allocate(childtemp % var) |
---|
84 | C |
---|
85 | childtemp % var % root_var => child % var % root_var |
---|
86 | C |
---|
87 | C Values of the grid variable |
---|
88 | childtemp % var % parray1 => tab |
---|
89 | C |
---|
90 | C Temporary results for the time interpolation before and after the space |
---|
91 | C interpolation |
---|
92 | childtemp % var % oldvalues2D => child % var % oldvalues2D |
---|
93 | C |
---|
94 | C Index indicating if a space interpolation is necessary |
---|
95 | childtemp % var % interpIndex => child % var % interpIndex |
---|
96 | childtemp % var % Interpolationshouldbemade = |
---|
97 | & child % var % Interpolationshouldbemade |
---|
98 | childtemp % var % list_interp => child % var% list_interp |
---|
99 | |
---|
100 | childtemp % var% lb = child % var % lb |
---|
101 | childtemp % var% ub = child % var % ub |
---|
102 | C |
---|
103 | C Call to the procedure for the calculations of the boundary conditions |
---|
104 | IF (present(procname)) THEN |
---|
105 | Call Agrif_CorrectVariable |
---|
106 | & (TypeInterp,parent,childtemp,deb,fin,pweight,weight,procname) |
---|
107 | ELSE |
---|
108 | Call Agrif_CorrectVariable |
---|
109 | & (TypeInterp,parent,childtemp,deb,fin,pweight,weight) |
---|
110 | ENDIF |
---|
111 | |
---|
112 | C |
---|
113 | child % var % oldvalues2D => childtemp % var % oldvalues2D |
---|
114 | child % var % list_interp => childtemp % var %list_interp |
---|
115 | C |
---|
116 | deallocate(childtemp % var) |
---|
117 | C |
---|
118 | C |
---|
119 | End Subroutine Agrif_Interp_bc_1D |
---|
120 | C |
---|
121 | |
---|
122 | C |
---|
123 | C |
---|
124 | C |
---|
125 | C ************************************************************************** |
---|
126 | CCC Subroutine Agrif_Interp_bc_2d |
---|
127 | C ************************************************************************** |
---|
128 | C |
---|
129 | Subroutine Agrif_Interp_bc_2d(TypeInterp,parent,child,tab,deb,fin, |
---|
130 | & weight,pweight,procname) |
---|
131 | C |
---|
132 | CCC Description: |
---|
133 | CCC Subroutine to calculate the boundary conditions on a fine grid for a 2D |
---|
134 | CCC grid variable. |
---|
135 | C |
---|
136 | C Declarations: |
---|
137 | C |
---|
138 | |
---|
139 | C |
---|
140 | C Arguments |
---|
141 | External :: procname |
---|
142 | Optional :: procname |
---|
143 | INTEGER,DIMENSION(6,6) :: TypeInterp ! TYPE of interpolation (linear, |
---|
144 | ! lagrange, spline, ... ) |
---|
145 | TYPE(AGRIF_PVariable) :: parent ! Variable on the parent grid |
---|
146 | TYPE(AGRIF_PVariable) :: child ! Variable on the child grid |
---|
147 | TYPE(AGRIF_PVariable) :: childtemp ! Temporary variable on the child |
---|
148 | ! grid |
---|
149 | INTEGER :: deb,fin ! Positions where interpolations are |
---|
150 | ! done on the fine grid |
---|
151 | REAL, DIMENSION( |
---|
152 | & child%var%lb(1):child%var%ub(1), |
---|
153 | & child%var%lb(2):child%var%ub(2) |
---|
154 | & ), Target :: tab ! Values of the grid variable |
---|
155 | LOGICAL :: pweight ! Indicates if weight is used for |
---|
156 | ! the temporal interpolation |
---|
157 | REAL :: weight ! Coefficient for the time |
---|
158 | ! interpolation |
---|
159 | C |
---|
160 | C |
---|
161 | C Definition of a temporary AGRIF_PVariable data TYPE representing the grid |
---|
162 | C variable. |
---|
163 | C |
---|
164 | allocate(childtemp % var) |
---|
165 | C |
---|
166 | childtemp % var % root_var => child % var % root_var |
---|
167 | C |
---|
168 | C Values of the grid variable |
---|
169 | childtemp % var % parray2 => tab |
---|
170 | C |
---|
171 | C Temporary results for the time interpolation before and after the space |
---|
172 | C interpolation |
---|
173 | childtemp % var % oldvalues2D => child % var % oldvalues2D |
---|
174 | C |
---|
175 | C Index indicating if a space interpolation is necessary |
---|
176 | childtemp % var % interpIndex => child % var % interpIndex |
---|
177 | childtemp % var % Interpolationshouldbemade = |
---|
178 | & child % var % Interpolationshouldbemade |
---|
179 | childtemp % var % list_interp => child % var% list_interp |
---|
180 | |
---|
181 | childtemp % var% lb = child % var % lb |
---|
182 | childtemp % var% ub = child % var % ub |
---|
183 | C |
---|
184 | C Call to the procedure for the calculations of the boundary conditions |
---|
185 | IF (present(procname)) THEN |
---|
186 | Call Agrif_CorrectVariable |
---|
187 | & (TypeInterp,parent,childtemp,deb,fin,pweight,weight,procname) |
---|
188 | ELSE |
---|
189 | Call Agrif_CorrectVariable |
---|
190 | & (TypeInterp,parent,childtemp,deb,fin,pweight,weight) |
---|
191 | ENDIF |
---|
192 | |
---|
193 | C |
---|
194 | child % var % oldvalues2D => childtemp % var % oldvalues2D |
---|
195 | child % var % list_interp => childtemp % var %list_interp |
---|
196 | C |
---|
197 | deallocate(childtemp % var) |
---|
198 | C |
---|
199 | C |
---|
200 | End Subroutine Agrif_Interp_bc_2D |
---|
201 | C |
---|
202 | C |
---|
203 | C |
---|
204 | C ************************************************************************** |
---|
205 | CCC Subroutine Agrif_Interp_bc_3d |
---|
206 | C ************************************************************************** |
---|
207 | C |
---|
208 | Subroutine Agrif_Interp_bc_3d(TypeInterp,parent,child,tab,deb,fin, |
---|
209 | & weight,pweight,procname) |
---|
210 | C |
---|
211 | CCC Description: |
---|
212 | CCC Subroutine to calculate the boundary conditions on a fine grid for a 3D |
---|
213 | CCC variable. |
---|
214 | C |
---|
215 | C Declarations: |
---|
216 | C |
---|
217 | |
---|
218 | C |
---|
219 | C Arguments |
---|
220 | External :: procname |
---|
221 | Optional :: procname |
---|
222 | INTEGER,DIMENSION(6,6) :: TypeInterp ! TYPE of interpolation (linear, |
---|
223 | ! lagrange, spline, ... ) |
---|
224 | TYPE(AGRIF_PVariable) :: parent ! Variable on the parent grid |
---|
225 | TYPE(AGRIF_PVariable) :: child ! Variable on the child grid |
---|
226 | TYPE(AGRIF_PVariable) :: childtemp ! Temporary variable on the child |
---|
227 | ! grid |
---|
228 | INTEGER :: deb,fin ! Positions where interpolations |
---|
229 | ! are done on the fine grid |
---|
230 | REAL, DIMENSION( |
---|
231 | & child%var%lb(1):child%var%ub(1), |
---|
232 | & child%var%lb(2):child%var%ub(2), |
---|
233 | & child%var%lb(3):child%var%ub(3) |
---|
234 | & ), Target :: tab ! Values of the grid variable |
---|
235 | LOGICAL :: pweight ! Indicates if weight is used for |
---|
236 | ! the temporal interpolation |
---|
237 | REAL :: weight ! Coefficient for the time |
---|
238 | ! interpolation |
---|
239 | C |
---|
240 | C |
---|
241 | C Definition of a temporary AGRIF_PVariable data TYPE representing the grid |
---|
242 | C variable. |
---|
243 | C |
---|
244 | allocate(childtemp % var) |
---|
245 | C |
---|
246 | childtemp % var % root_var => child % var % root_var |
---|
247 | C |
---|
248 | C Values of the grid variable |
---|
249 | childtemp % var % parray3 => tab |
---|
250 | C |
---|
251 | C Temporary results for the time interpolation before and after the space |
---|
252 | C interpolation |
---|
253 | childtemp % var % oldvalues2D => child % var % oldvalues2D |
---|
254 | C |
---|
255 | C Index indicating if a space interpolation is necessary |
---|
256 | childtemp % var % interpIndex => child % var % interpIndex |
---|
257 | childtemp % var % Interpolationshouldbemade = |
---|
258 | & child % var % Interpolationshouldbemade |
---|
259 | childtemp % var % list_interp => child % var% list_interp |
---|
260 | |
---|
261 | childtemp % var% lb = child % var % lb |
---|
262 | childtemp % var% ub = child % var % ub |
---|
263 | C |
---|
264 | C Call to the procedure for the calculations of the boundary conditions |
---|
265 | IF (present(procname)) THEN |
---|
266 | Call Agrif_CorrectVariable |
---|
267 | & (TypeInterp,parent,childtemp,deb,fin,pweight,weight,procname) |
---|
268 | ELSE |
---|
269 | Call Agrif_CorrectVariable |
---|
270 | & (TypeInterp,parent,childtemp,deb,fin,pweight,weight) |
---|
271 | ENDIF |
---|
272 | C |
---|
273 | child % var % oldvalues2D => childtemp % var % oldvalues2D |
---|
274 | child % var % list_interp => childtemp % var %list_interp |
---|
275 | C |
---|
276 | deallocate(childtemp % var) |
---|
277 | C |
---|
278 | C |
---|
279 | End Subroutine Agrif_Interp_bc_3D |
---|
280 | C |
---|
281 | C |
---|
282 | C |
---|
283 | C |
---|
284 | C ************************************************************************** |
---|
285 | CCC Subroutine Agrif_Interp_bc_4d |
---|
286 | C ************************************************************************** |
---|
287 | C |
---|
288 | Subroutine Agrif_Interp_bc_4d(TypeInterp,parent,child,tab,deb,fin, |
---|
289 | & weight,pweight,procname) |
---|
290 | C |
---|
291 | CCC Description: |
---|
292 | CCC Subroutine to calculate the boundary conditions on a fine grid for a 4D |
---|
293 | CCC grid variable. |
---|
294 | C |
---|
295 | C Declarations: |
---|
296 | C |
---|
297 | |
---|
298 | C |
---|
299 | C Arguments |
---|
300 | External :: procname |
---|
301 | Optional :: procname |
---|
302 | INTEGER,DIMENSION(6,6) :: TypeInterp ! TYPE of interpolation (linear, |
---|
303 | ! lagrange, spline, ... ) |
---|
304 | TYPE(AGRIF_PVariable) :: parent ! Variable on the parent grid |
---|
305 | TYPE(AGRIF_PVariable) :: child ! Variable on the child grid |
---|
306 | TYPE(AGRIF_PVariable) :: childtemp ! Temporary varaiable on the child |
---|
307 | ! grid |
---|
308 | INTEGER :: deb,fin ! Positions where interpolations |
---|
309 | ! are done on the fine grid |
---|
310 | REAL, DIMENSION( |
---|
311 | & child%var%lb(1):child%var%ub(1), |
---|
312 | & child%var%lb(2):child%var%ub(2), |
---|
313 | & child%var%lb(3):child%var%ub(3), |
---|
314 | & child%var%lb(4):child%var%ub(4) |
---|
315 | & ), Target :: tab ! Values of the grid variable |
---|
316 | LOGICAL :: pweight ! Indicates if weight is used for |
---|
317 | ! the temporal interpolation |
---|
318 | REAL :: weight ! Coefficient for the time |
---|
319 | ! interpolation |
---|
320 | C |
---|
321 | C |
---|
322 | C Definition of a temporary AGRIF_PVariable data TYPE representing the grid |
---|
323 | C variable. |
---|
324 | C |
---|
325 | allocate(childtemp % var) |
---|
326 | C |
---|
327 | childtemp % var % root_var => child % var % root_var |
---|
328 | C |
---|
329 | C Values of the grid variable |
---|
330 | childtemp % var % parray4 => tab |
---|
331 | C |
---|
332 | C Temporary results for the time interpolation before and after the space |
---|
333 | C interpolation |
---|
334 | childtemp % var % oldvalues2D => child % var % oldvalues2D |
---|
335 | C |
---|
336 | C Index indicating if a space interpolation is necessary |
---|
337 | childtemp % var % interpIndex => child % var % interpIndex |
---|
338 | childtemp % var % Interpolationshouldbemade = |
---|
339 | & child % var % Interpolationshouldbemade |
---|
340 | childtemp % var % list_interp => child % var% list_interp |
---|
341 | |
---|
342 | childtemp % var% lb = child % var % lb |
---|
343 | childtemp % var% ub = child % var % ub |
---|
344 | C |
---|
345 | C Call to the procedure for the calculations of the boundary conditions |
---|
346 | IF (present(procname)) THEN |
---|
347 | Call Agrif_CorrectVariable |
---|
348 | & (TypeInterp,parent,childtemp,deb,fin,pweight,weight,procname) |
---|
349 | ELSE |
---|
350 | Call Agrif_CorrectVariable |
---|
351 | & (TypeInterp,parent,childtemp,deb,fin,pweight,weight) |
---|
352 | ENDIF |
---|
353 | C |
---|
354 | child % var % oldvalues2D => childtemp % var % oldvalues2D |
---|
355 | child % var % list_interp => childtemp % var %list_interp |
---|
356 | C |
---|
357 | deallocate(childtemp % var) |
---|
358 | C |
---|
359 | C |
---|
360 | End Subroutine Agrif_Interp_bc_4D |
---|
361 | C |
---|
362 | C |
---|
363 | C |
---|
364 | C ************************************************************************** |
---|
365 | CCC Subroutine Agrif_Interp_bc_5d |
---|
366 | C ************************************************************************** |
---|
367 | C |
---|
368 | Subroutine Agrif_Interp_bc_5d(TypeInterp,parent,child,tab,deb,fin, |
---|
369 | & weight,pweight,procname) |
---|
370 | C |
---|
371 | CCC Description: |
---|
372 | CCC Subroutine to calculate the boundary conditions on a fine grid for a 5D |
---|
373 | CCC grid variable. |
---|
374 | C |
---|
375 | C Declarations: |
---|
376 | C |
---|
377 | |
---|
378 | C |
---|
379 | C Arguments |
---|
380 | External :: procname |
---|
381 | Optional :: procname |
---|
382 | INTEGER,DIMENSION(6,6) :: TypeInterp ! TYPE of interpolation (linear, |
---|
383 | ! lagrange, spline, ... ) |
---|
384 | TYPE(AGRIF_PVariable) :: parent ! Variable on the parent grid |
---|
385 | TYPE(AGRIF_PVariable) :: child ! Variable on the child grid |
---|
386 | TYPE(AGRIF_PVariable) :: childtemp ! Temporary varaiable on the child |
---|
387 | ! grid |
---|
388 | INTEGER :: deb,fin ! Positions where interpolations |
---|
389 | ! are done on the fine grid |
---|
390 | REAL, DIMENSION( |
---|
391 | & child%var%lb(1):child%var%ub(1), |
---|
392 | & child%var%lb(2):child%var%ub(2), |
---|
393 | & child%var%lb(3):child%var%ub(3), |
---|
394 | & child%var%lb(4):child%var%ub(4), |
---|
395 | & child%var%lb(5):child%var%ub(5) |
---|
396 | & ), Target :: tab ! Values of the grid variable |
---|
397 | LOGICAL :: pweight ! Indicates if weight is used for |
---|
398 | ! the temporal interpolation |
---|
399 | REAL :: weight ! Coefficient for the time |
---|
400 | ! interpolation |
---|
401 | C |
---|
402 | C |
---|
403 | C Definition of a temporary AGRIF_PVariable data TYPE representing the grid |
---|
404 | C variable. |
---|
405 | C |
---|
406 | allocate(childtemp % var) |
---|
407 | C |
---|
408 | childtemp % var % root_var => child % var % root_var |
---|
409 | C |
---|
410 | C Values of the grid variable |
---|
411 | childtemp % var % parray5 => tab |
---|
412 | C |
---|
413 | C Temporary results for the time interpolation before and after the space |
---|
414 | C interpolation |
---|
415 | childtemp % var % oldvalues2D => child % var % oldvalues2D |
---|
416 | C |
---|
417 | C Index indicating if a space interpolation is necessary |
---|
418 | childtemp % var % interpIndex => child % var % interpIndex |
---|
419 | childtemp % var % Interpolationshouldbemade = |
---|
420 | & child % var % Interpolationshouldbemade |
---|
421 | childtemp % var % list_interp => child % var% list_interp |
---|
422 | |
---|
423 | childtemp % var% lb = child % var % lb |
---|
424 | childtemp % var% ub = child % var % ub |
---|
425 | |
---|
426 | C |
---|
427 | C Call to the procedure for the calculations of the boundary conditions |
---|
428 | IF (present(procname)) THEN |
---|
429 | Call Agrif_CorrectVariable |
---|
430 | & (TypeInterp,parent,childtemp,deb,fin,pweight,weight,procname) |
---|
431 | ELSE |
---|
432 | Call Agrif_CorrectVariable |
---|
433 | & (TypeInterp,parent,childtemp,deb,fin,pweight,weight) |
---|
434 | ENDIF |
---|
435 | C |
---|
436 | child % var % oldvalues2D => childtemp % var % oldvalues2D |
---|
437 | child % var % list_interp => childtemp % var %list_interp |
---|
438 | C |
---|
439 | deallocate(childtemp % var) |
---|
440 | C |
---|
441 | C |
---|
442 | End Subroutine Agrif_Interp_bc_5D |
---|
443 | C |
---|
444 | C |
---|
445 | C |
---|
446 | C |
---|
447 | C ************************************************************************** |
---|
448 | CCC Subroutine Agrif_Interp_bc_6d |
---|
449 | C ************************************************************************** |
---|
450 | C |
---|
451 | Subroutine Agrif_Interp_bc_6d(TypeInterp,parent,child,tab,deb,fin, |
---|
452 | & weight,pweight) |
---|
453 | C |
---|
454 | CCC Description: |
---|
455 | CCC Subroutine to calculate the boundary conditions on a fine grid for a 6D |
---|
456 | CCC grid variable. |
---|
457 | C |
---|
458 | C Declarations: |
---|
459 | C |
---|
460 | |
---|
461 | C |
---|
462 | C Arguments |
---|
463 | INTEGER,DIMENSION(6,6) :: TypeInterp ! TYPE of interpolation (linear, |
---|
464 | ! lagrange, spline, ... ) |
---|
465 | TYPE(AGRIF_PVariable) :: parent ! Variable on the parent grid |
---|
466 | TYPE(AGRIF_PVariable) :: child ! Variable on the child grid |
---|
467 | TYPE(AGRIF_PVariable) :: childtemp ! Temporary varaiable on the child |
---|
468 | ! grid |
---|
469 | INTEGER :: deb,fin ! Positions where interpolations |
---|
470 | ! are done on the fine grid |
---|
471 | REAL, DIMENSION( |
---|
472 | & child%var%lb(1):child%var%ub(1), |
---|
473 | & child%var%lb(2):child%var%ub(2), |
---|
474 | & child%var%lb(3):child%var%ub(3), |
---|
475 | & child%var%lb(4):child%var%ub(4), |
---|
476 | & child%var%lb(5):child%var%ub(5), |
---|
477 | & child%var%lb(6):child%var%ub(6) |
---|
478 | & ), Target :: tab ! Values of the grid variable |
---|
479 | LOGICAL :: pweight ! Indicates if weight is used for |
---|
480 | ! the temporal interpolation |
---|
481 | REAL :: weight ! Coefficient for the time |
---|
482 | ! interpolation |
---|
483 | C |
---|
484 | C |
---|
485 | C Definition of a temporary AGRIF_PVariable data TYPE representing the grid |
---|
486 | C variable. |
---|
487 | C |
---|
488 | allocate(childtemp % var) |
---|
489 | C |
---|
490 | childtemp % var % root_var => child % var % root_var |
---|
491 | C |
---|
492 | C Values of the grid variable |
---|
493 | childtemp % var % parray6 => tab |
---|
494 | C |
---|
495 | C Temporary results for the time interpolation before and after the space |
---|
496 | C interpolation |
---|
497 | childtemp % var % oldvalues2D => child % var % oldvalues2D |
---|
498 | C |
---|
499 | C Index indicating if a space interpolation is necessary |
---|
500 | childtemp % var % interpIndex => child % var % interpIndex |
---|
501 | childtemp % var % Interpolationshouldbemade = |
---|
502 | & child % var % Interpolationshouldbemade |
---|
503 | childtemp % var % list_interp => child % var% list_interp |
---|
504 | |
---|
505 | childtemp % var% lb = child % var % lb |
---|
506 | childtemp % var% ub = child % var % ub |
---|
507 | C |
---|
508 | C Call to the procedure for the calculations of the boundary conditions |
---|
509 | Call Agrif_CorrectVariable |
---|
510 | & (TypeInterp,parent,childtemp,deb,fin,pweight,weight) |
---|
511 | C |
---|
512 | child % var % oldvalues2D => childtemp % var % oldvalues2D |
---|
513 | child % var % list_interp => childtemp % var %list_interp |
---|
514 | C |
---|
515 | deallocate(childtemp % var) |
---|
516 | C |
---|
517 | C |
---|
518 | End Subroutine Agrif_Interp_bc_6D |
---|
519 | C |
---|
520 | C |
---|
521 | C ************************************************************************** |
---|
522 | CCC Subroutine Agrif_CorrectVariable |
---|
523 | C ************************************************************************** |
---|
524 | C |
---|
525 | Subroutine AGRIF_CorrectVariable(TypeInterp,parent,child,deb,fin, |
---|
526 | & pweight,weight,procname) |
---|
527 | C |
---|
528 | CCC Description: |
---|
529 | CCC Subroutine to calculate the boundary conditions on a fine grid. |
---|
530 | C |
---|
531 | C Declarations: |
---|
532 | C |
---|
533 | |
---|
534 | C |
---|
535 | C Arguments |
---|
536 | External :: procname |
---|
537 | Optional :: procname |
---|
538 | TYPE(AGRIF_PVariable) :: parent ! Variable on the parent grid |
---|
539 | TYPE(AGRIF_PVariable) :: child ! Variable on the child grid |
---|
540 | INTEGER,DIMENSION(6,6) :: TypeInterp ! TYPE of interpolation |
---|
541 | ! (linear,lagrange,...) |
---|
542 | INTEGER :: deb,fin ! Positions where boundary |
---|
543 | ! conditions are calculated |
---|
544 | LOGICAL :: pweight ! Indicates if weight is used |
---|
545 | ! for the time interpolation |
---|
546 | REAL :: weight ! Coefficient for the time |
---|
547 | ! interpolation |
---|
548 | C |
---|
549 | C Local scalars |
---|
550 | TYPE(Agrif_Grid) , Pointer :: Agrif_Child_Gr,Agrif_Parent_Gr |
---|
551 | TYPE(AGRIF_Variable), Pointer :: root ! Variable on the root grid |
---|
552 | INTEGER :: nbdim ! Number of dimensions of |
---|
553 | ! the grid variable |
---|
554 | INTEGER :: n |
---|
555 | INTEGER,DIMENSION(6) :: pttab_child ! Index of the first point |
---|
556 | ! inside the domain for |
---|
557 | ! the child grid variable |
---|
558 | INTEGER,DIMENSION(6) :: pttab_parent ! Index of the first point |
---|
559 | ! inside the domain for |
---|
560 | ! the parent grid |
---|
561 | ! variable |
---|
562 | INTEGER,DIMENSION(6) :: nbtab_Child ! Number of the cells |
---|
563 | INTEGER,DIMENSION(6) :: posvartab_Child ! Position of the |
---|
564 | ! variable on the cell |
---|
565 | INTEGER,DIMENSION(6) :: loctab_Child ! Indicates if the child |
---|
566 | ! grid has a common |
---|
567 | ! border with the root |
---|
568 | ! grid |
---|
569 | REAL, DIMENSION(6) :: s_child,s_parent ! Positions of the |
---|
570 | ! parent and child grids |
---|
571 | REAL, DIMENSION(6) :: ds_child,ds_parent ! Space steps of the |
---|
572 | ! parent and child grids |
---|
573 | C |
---|
574 | C |
---|
575 | loctab_child(:) = 0 |
---|
576 | C |
---|
577 | Agrif_Child_Gr => Agrif_Curgrid |
---|
578 | Agrif_Parent_Gr => Agrif_Curgrid % parent |
---|
579 | root => child % var % root_var |
---|
580 | nbdim = root % nbdim |
---|
581 | C |
---|
582 | do n = 1,nbdim |
---|
583 | posvartab_child(n) = root % posvar(n) |
---|
584 | enddo |
---|
585 | C |
---|
586 | C |
---|
587 | do n = 1,nbdim |
---|
588 | C |
---|
589 | Select case(root % interptab(n)) |
---|
590 | C |
---|
591 | case('x') ! x DIMENSION |
---|
592 | C |
---|
593 | nbtab_Child(n) = Agrif_Child_Gr % nb(1) |
---|
594 | pttab_Child(n) = root % point(1) |
---|
595 | pttab_Parent(n) = root % point(1) |
---|
596 | s_Child(n) = Agrif_Child_Gr % Agrif_x(1) |
---|
597 | s_Parent(n) = Agrif_Parent_Gr % Agrif_x(1) |
---|
598 | ds_Child(n) = Agrif_Child_Gr % Agrif_d(1) |
---|
599 | ds_Parent(n) = Agrif_Parent_Gr % Agrif_d(1) |
---|
600 | if (root % posvar(n) == 2) then |
---|
601 | s_Child(n) = s_Child(n) + ds_Child(n)/2. |
---|
602 | s_Parent(n) = s_Parent(n) + ds_Parent(n)/2. |
---|
603 | endif |
---|
604 | C |
---|
605 | if (Agrif_CURGRID % NearRootBorder(1)) |
---|
606 | & loctab_child(n) = -1 |
---|
607 | if (Agrif_CURGRID % DistantRootBorder(1)) |
---|
608 | & loctab_child(n) = -2 |
---|
609 | if ((Agrif_CURGRID % NearRootBorder(1)) .AND. |
---|
610 | & (Agrif_CURGRID % DistantRootBorder(1))) |
---|
611 | & loctab_child(n) = -3 |
---|
612 | C |
---|
613 | case('y') ! y DIMENSION |
---|
614 | C |
---|
615 | nbtab_Child(n) = Agrif_Child_Gr % nb(2) |
---|
616 | pttab_Child(n) = root % point(2) |
---|
617 | pttab_Parent(n) = root % point(2) |
---|
618 | s_Child(n) = Agrif_Child_Gr % Agrif_x(2) |
---|
619 | s_Parent(n) = Agrif_Parent_Gr % Agrif_x(2) |
---|
620 | ds_Child(n) = Agrif_Child_Gr % Agrif_d(2) |
---|
621 | ds_Parent(n) = Agrif_Parent_Gr % Agrif_d(2) |
---|
622 | if (root % posvar(n) == 2) then |
---|
623 | s_Child(n) = s_Child(n) + ds_Child(n)/2. |
---|
624 | s_Parent(n) = s_Parent(n) + ds_Parent(n)/2. |
---|
625 | endif |
---|
626 | C |
---|
627 | if (Agrif_CURGRID % NearRootBorder(2)) |
---|
628 | & loctab_child(n) = -1 |
---|
629 | if (Agrif_CURGRID % DistantRootBorder(2)) |
---|
630 | & loctab_child(n) = -2 |
---|
631 | if ((Agrif_CURGRID % NearRootBorder(2)) .AND. |
---|
632 | & (Agrif_CURGRID % DistantRootBorder(2))) |
---|
633 | & loctab_child(n) = -3 |
---|
634 | C |
---|
635 | case('z') ! z DIMENSION |
---|
636 | C |
---|
637 | nbtab_Child(n) = Agrif_Child_Gr % nb(3) |
---|
638 | pttab_Child(n) = root % point(3) |
---|
639 | pttab_Parent(n) = root % point(3) |
---|
640 | s_Child(n) = Agrif_Child_Gr % Agrif_x(3) |
---|
641 | s_Parent(n) = Agrif_Parent_Gr % Agrif_x(3) |
---|
642 | ds_Child(n) = Agrif_Child_Gr % Agrif_d(3) |
---|
643 | ds_Parent(n) = Agrif_Parent_Gr % Agrif_d(3) |
---|
644 | if (root % posvar(n) == 2) then |
---|
645 | s_Child(n) = s_Child(n) + ds_Child(n)/2. |
---|
646 | s_Parent(n) = s_Parent(n) + ds_Parent(n)/2. |
---|
647 | endif |
---|
648 | C |
---|
649 | if (Agrif_CURGRID % NearRootBorder(3)) |
---|
650 | & loctab_child(n) = -1 |
---|
651 | if (Agrif_CURGRID % DistantRootBorder(3)) |
---|
652 | & loctab_child(n) = -2 |
---|
653 | if ((Agrif_CURGRID % NearRootBorder(3)) .AND. |
---|
654 | & (Agrif_CURGRID % DistantRootBorder(3))) |
---|
655 | & loctab_child(n) = -3 |
---|
656 | C |
---|
657 | case('N') ! No space DIMENSION |
---|
658 | C |
---|
659 | nbtab_Child(n) = child % var % ub(n) - child % var % lb(n) |
---|
660 | pttab_Child(n) = child % var % lb(n) |
---|
661 | C |
---|
662 | C No interpolation but only a copy of the values of the grid variable |
---|
663 | C |
---|
664 | posvartab_child(n) = 1 |
---|
665 | pttab_Parent(n)= pttab_Child(n) |
---|
666 | s_Child(n)=0. |
---|
667 | s_Parent(n)=0. |
---|
668 | ds_Child(n)=1. |
---|
669 | ds_Parent(n)=1. |
---|
670 | loctab_child(n) = -3 |
---|
671 | C |
---|
672 | End select |
---|
673 | C |
---|
674 | enddo |
---|
675 | C |
---|
676 | IF (present(procname)) THEN |
---|
677 | Call AGRIF_CorrectnD |
---|
678 | & (TypeInterp,parent,child,deb,fin,pweight,weight, |
---|
679 | & pttab_Child(1:nbdim),pttab_Parent(1:nbdim), |
---|
680 | & nbtab_Child(1:nbdim),posvartab_Child(1:nbdim), |
---|
681 | & loctab_Child(1:nbdim), |
---|
682 | & s_Child(1:nbdim),s_Parent(1:nbdim), |
---|
683 | & ds_Child(1:nbdim),ds_Parent(1:nbdim),nbdim,procname) |
---|
684 | ELSE |
---|
685 | Call AGRIF_CorrectnD |
---|
686 | & (TypeInterp,parent,child,deb,fin,pweight,weight, |
---|
687 | & pttab_Child(1:nbdim),pttab_Parent(1:nbdim), |
---|
688 | & nbtab_Child(1:nbdim),posvartab_Child(1:nbdim), |
---|
689 | & loctab_Child(1:nbdim), |
---|
690 | & s_Child(1:nbdim),s_Parent(1:nbdim), |
---|
691 | & ds_Child(1:nbdim),ds_Parent(1:nbdim),nbdim) |
---|
692 | ENDIF |
---|
693 | C |
---|
694 | C |
---|
695 | End subroutine AGRIF_CorrectVariable |
---|
696 | C |
---|
697 | C ************************************************************************** |
---|
698 | CCC Subroutine Agrif_Correctnd |
---|
699 | C ************************************************************************** |
---|
700 | C |
---|
701 | Subroutine AGRIF_Correctnd(TypeInterp,parent,child,deb,fin, |
---|
702 | & pweight,weight, |
---|
703 | & pttab_child,pttab_Parent, |
---|
704 | & nbtab_Child,posvartab_Child, |
---|
705 | & loctab_Child, |
---|
706 | & s_Child,s_Parent, |
---|
707 | & ds_Child,ds_Parent,nbdim,procname) |
---|
708 | C |
---|
709 | CCC Description: |
---|
710 | CCC Subroutine to calculate the boundary conditions for a nD grid variable on |
---|
711 | CCC a fine grid by using a space and time interpolations; it is called by the |
---|
712 | CCC Agrif_CorrectVariable procedure. |
---|
713 | C |
---|
714 | C |
---|
715 | C Declarations: |
---|
716 | C |
---|
717 | |
---|
718 | C |
---|
719 | #ifdef key_mpp_mpi |
---|
720 | C |
---|
721 | INCLUDE 'mpif.h' |
---|
722 | C |
---|
723 | #endif |
---|
724 | C |
---|
725 | C Arguments |
---|
726 | External :: procname |
---|
727 | Optional :: procname |
---|
728 | INTEGER,DIMENSION(6,6) :: TypeInterp ! TYPE of interpolation (linear, |
---|
729 | ! spline,...) |
---|
730 | TYPE(AGRIF_PVariable) :: parent ! Variable on the parent grid |
---|
731 | TYPE(AGRIF_PVariable) :: child ! Variable on the child grid |
---|
732 | INTEGER :: deb,fin ! Positions where interpolations |
---|
733 | ! are done |
---|
734 | LOGICAL :: pweight ! Indicates if weight is used for |
---|
735 | ! the temporal interpolation |
---|
736 | REAL :: weight ! Coefficient for the temporal |
---|
737 | ! interpolation |
---|
738 | INTEGER :: nbdim ! Number of dimensions of the grid |
---|
739 | ! variable |
---|
740 | INTEGER,DIMENSION(nbdim) :: pttab_child ! Index of the first point inside |
---|
741 | ! the domain for the parent |
---|
742 | ! grid variable |
---|
743 | INTEGER,DIMENSION(nbdim) :: pttab_Parent ! Index of the first point |
---|
744 | ! inside the domain for the |
---|
745 | ! child grid variable |
---|
746 | INTEGER,DIMENSION(nbdim) :: nbtab_Child ! Number of cells of the child |
---|
747 | ! grid |
---|
748 | INTEGER,DIMENSION(nbdim) :: posvartab_Child ! Position of the grid |
---|
749 | ! variable (1 or 2) |
---|
750 | INTEGER,DIMENSION(nbdim) :: loctab_Child ! Indicates if the child |
---|
751 | ! grid has a common border with |
---|
752 | ! the root grid |
---|
753 | REAL ,DIMENSION(nbdim) :: s_Child,s_Parent ! Positions of the parent |
---|
754 | ! and child grids |
---|
755 | REAL ,DIMENSION(nbdim) :: ds_Child,ds_Parent ! Space steps of the |
---|
756 | ! parent and child grids |
---|
757 | C |
---|
758 | C Local variables |
---|
759 | TYPE(AGRIF_PVariable) :: restore ! Variable on the parent |
---|
760 | INTEGER,DIMENSION(nbdim,2) :: lubglob |
---|
761 | INTEGER :: i |
---|
762 | INTEGER :: kindex ! Index used for safeguard |
---|
763 | ! and time interpolation |
---|
764 | INTEGER,DIMENSION(nbdim,2,2) :: indtab ! Arrays indicating the limits |
---|
765 | ! of the child |
---|
766 | INTEGER,DIMENSION(nbdim,2,2) :: indtruetab ! grid variable where |
---|
767 | ! boundary conditions are |
---|
768 | INTEGER,DIMENSION(nbdim,2,2,nbdim) :: ptres,ptres2 ! calculated |
---|
769 | INTEGER :: nb,ndir,n,sizetab(1) |
---|
770 | REAL, DIMENSION(:), Allocatable :: tab ! Array used for the interpolation |
---|
771 | REAL :: c1t,c2t ! Coefficients for the time interpolation |
---|
772 | ! (c2t=1-c1t) |
---|
773 | C |
---|
774 | #ifdef key_mpp_mpi |
---|
775 | C |
---|
776 | INTEGER,DIMENSION(nbdim) :: lower,upper |
---|
777 | INTEGER,DIMENSION(nbdim) :: ltab,utab |
---|
778 | INTEGER,DIMENSION(nbdim) :: lb,ub |
---|
779 | INTEGER,DIMENSION(nbdim,2) :: iminmaxg |
---|
780 | INTEGER :: code |
---|
781 | C |
---|
782 | #endif |
---|
783 | C |
---|
784 | C |
---|
785 | indtab(1:nbdim,2,1) = pttab_child(1:nbdim) + nbtab_child(1:nbdim) |
---|
786 | & + deb |
---|
787 | indtab(1:nbdim,2,2) = indtab(1:nbdim,2,1) + ( fin - deb ) |
---|
788 | |
---|
789 | indtab(1:nbdim,1,1) = pttab_child(1:nbdim) - fin |
---|
790 | indtab(1:nbdim,1,2) = pttab_child(1:nbdim) - deb |
---|
791 | |
---|
792 | WHERE (posvartab_child(1:nbdim) == 2) |
---|
793 | indtab(1:nbdim,1,1) = indtab(1:nbdim,1,1) - 1 |
---|
794 | indtab(1:nbdim,1,2) = indtab(1:nbdim,1,2) - 1 |
---|
795 | END WHERE |
---|
796 | |
---|
797 | |
---|
798 | #if !defined key_mpp_mpi |
---|
799 | Call Agrif_nbdim_Get_bound_dimension(child%var,lubglob(:,1), |
---|
800 | & lubglob(:,2),nbdim) |
---|
801 | C |
---|
802 | #else |
---|
803 | C |
---|
804 | Call Agrif_nbdim_Get_bound_dimension(child%var,lb,ub,nbdim) |
---|
805 | |
---|
806 | DO i = 1,nbdim |
---|
807 | C |
---|
808 | Call Agrif_Invloc(lb(i),Agrif_Procrank,i,iminmaxg(i,1)) |
---|
809 | Call Agrif_Invloc(ub(i),Agrif_Procrank,i,iminmaxg(i,2)) |
---|
810 | C |
---|
811 | ENDDO |
---|
812 | C |
---|
813 | iminmaxg(1:nbdim,2) = - iminmaxg(1:nbdim,2) |
---|
814 | |
---|
815 | CALL MPI_ALLREDUCE(iminmaxg,lubglob,2*nbdim,MPI_INTEGER,MPI_MIN, |
---|
816 | & MPI_COMM_AGRIF,code) |
---|
817 | |
---|
818 | lubglob(1:nbdim,2) = - lubglob(1:nbdim,2) |
---|
819 | C |
---|
820 | #endif |
---|
821 | C |
---|
822 | indtruetab(1:nbdim,1,1) = max(indtab(1:nbdim,1,1), |
---|
823 | & lubglob(1:nbdim,1)) |
---|
824 | indtruetab(1:nbdim,1,2) = max(indtab(1:nbdim,1,2), |
---|
825 | & lubglob(1:nbdim,1)) |
---|
826 | indtruetab(1:nbdim,2,1) = min(indtab(1:nbdim,2,1), |
---|
827 | & lubglob(1:nbdim,2)) |
---|
828 | indtruetab(1:nbdim,2,2) = min(indtab(1:nbdim,2,2), |
---|
829 | & lubglob(1:nbdim,2)) |
---|
830 | |
---|
831 | |
---|
832 | C |
---|
833 | C |
---|
834 | do nb = 1,nbdim |
---|
835 | C |
---|
836 | do ndir = 1,2 |
---|
837 | C |
---|
838 | if (loctab_child(nb) /= (-ndir) |
---|
839 | & .AND. loctab_child(nb) /= -3) then |
---|
840 | C |
---|
841 | do n = 1,2 |
---|
842 | C |
---|
843 | ptres(nb,n,ndir,nb) = indtruetab(nb,ndir,n) |
---|
844 | C |
---|
845 | enddo |
---|
846 | C |
---|
847 | do i = 1,nbdim |
---|
848 | C |
---|
849 | if (i .NE. nb) then |
---|
850 | C |
---|
851 | if (loctab_child(i) == -1 |
---|
852 | & .OR. loctab_child(i) == -3) then |
---|
853 | C |
---|
854 | ptres(i,1,ndir,nb) = pttab_child(i) |
---|
855 | C |
---|
856 | else |
---|
857 | C |
---|
858 | ptres(i,1,ndir,nb) = indtruetab(i,1,1) |
---|
859 | C |
---|
860 | endif |
---|
861 | C |
---|
862 | if (loctab_child(i) == -2 |
---|
863 | & .OR. loctab_child(i) == -3) then |
---|
864 | C |
---|
865 | if (posvartab_child(i) == 1) then |
---|
866 | C |
---|
867 | ptres(i,2,ndir,nb) = pttab_child(i) |
---|
868 | & + nbtab_child(i) |
---|
869 | C |
---|
870 | else |
---|
871 | C |
---|
872 | ptres(i,2,ndir,nb) = pttab_child(i) |
---|
873 | & + nbtab_child(i) - 1 |
---|
874 | C |
---|
875 | endif |
---|
876 | C |
---|
877 | else |
---|
878 | C |
---|
879 | ptres(i,2,ndir,nb) = indtruetab(i,2,2) |
---|
880 | C |
---|
881 | endif |
---|
882 | C |
---|
883 | endif |
---|
884 | C |
---|
885 | enddo |
---|
886 | |
---|
887 | C |
---|
888 | #if defined key_mpp_mpi |
---|
889 | Call Agrif_nbdim_Get_bound_dimension |
---|
890 | & (child%var,lower,upper,nbdim) |
---|
891 | |
---|
892 | do i = 1,nbdim |
---|
893 | C |
---|
894 | Call GetLocalBoundaries(ptres(i,1,ndir,nb), |
---|
895 | & ptres(i,2,ndir,nb),i, |
---|
896 | & lower(i),upper(i), |
---|
897 | & ltab(i),utab(i)) |
---|
898 | ptres2(i,1,ndir,nb) = max(ltab(i),lower(i)) |
---|
899 | ptres2(i,2,ndir,nb) = min(utab(i),upper(i)) |
---|
900 | if ((i == nb) .AND. (ndir == 1)) then |
---|
901 | ptres2(i,2,ndir,nb) = max(utab(i),lower(i)) |
---|
902 | elseif ((i == nb) .AND. (ndir == 2)) then |
---|
903 | ptres2(i,1,ndir,nb) = min(ltab(i),upper(i)) |
---|
904 | endif |
---|
905 | C |
---|
906 | enddo |
---|
907 | #else |
---|
908 | ptres2(:,:,ndir,nb) = ptres(:,:,ndir,nb) |
---|
909 | #endif |
---|
910 | |
---|
911 | endif |
---|
912 | |
---|
913 | enddo |
---|
914 | enddo |
---|
915 | C |
---|
916 | if (child % var % interpIndex |
---|
917 | & /= Agrif_Curgrid % parent % ngridstep .OR. |
---|
918 | & child%var%Interpolationshouldbemade ) then |
---|
919 | C |
---|
920 | C Space interpolation |
---|
921 | C |
---|
922 | kindex = 1 |
---|
923 | |
---|
924 | |
---|
925 | C |
---|
926 | do nb = 1,nbdim |
---|
927 | C |
---|
928 | do ndir = 1,2 |
---|
929 | C |
---|
930 | if (loctab_child(nb) /= (-ndir) |
---|
931 | & .AND. loctab_child(nb) /= -3) then |
---|
932 | C |
---|
933 | IF (present(procname)) THEN |
---|
934 | Call Agrif_InterpnD |
---|
935 | & (TYPEInterp(nb,:),parent,child, |
---|
936 | & ptres(1:nbdim,1,ndir,nb),ptres(1:nbdim,2,ndir,nb), |
---|
937 | & pttab_child(1:nbdim),pttab_Parent(1:nbdim), |
---|
938 | & s_Child(1:nbdim),s_Parent(1:nbdim), |
---|
939 | & ds_Child(1:nbdim),ds_Parent(1:nbdim), |
---|
940 | & restore,.FALSE.,nbdim,procname) |
---|
941 | ELSE |
---|
942 | Call Agrif_InterpnD |
---|
943 | & (TYPEInterp(nb,:),parent,child, |
---|
944 | & ptres(1:nbdim,1,ndir,nb),ptres(1:nbdim,2,ndir,nb), |
---|
945 | & pttab_child(1:nbdim),pttab_Parent(1:nbdim), |
---|
946 | & s_Child(1:nbdim),s_Parent(1:nbdim), |
---|
947 | & ds_Child(1:nbdim),ds_Parent(1:nbdim), |
---|
948 | & restore,.FALSE.,nbdim) |
---|
949 | ENDIF |
---|
950 | |
---|
951 | IF (.NOT. child%var%interpolationshouldbemade) THEN |
---|
952 | C |
---|
953 | C Safeguard of the values of the grid variable (at times n and n+1 |
---|
954 | C on the parent grid) |
---|
955 | C |
---|
956 | sizetab(1) = 1 |
---|
957 | C |
---|
958 | do i = 1,nbdim |
---|
959 | C |
---|
960 | sizetab(1) = sizetab(1) |
---|
961 | & * (ptres2(i,2,ndir,nb)-ptres2(i,1,ndir,nb)+1) |
---|
962 | C |
---|
963 | enddo |
---|
964 | |
---|
965 | Call saveAfterInterp(child, |
---|
966 | & ptres2(:,:,ndir,nb),kindex,sizetab(1),nbdim) |
---|
967 | C |
---|
968 | ENDIF |
---|
969 | C |
---|
970 | endif |
---|
971 | C |
---|
972 | enddo |
---|
973 | C |
---|
974 | enddo |
---|
975 | C |
---|
976 | C |
---|
977 | child % var % interpIndex = Agrif_Curgrid % parent % ngridstep |
---|
978 | C |
---|
979 | C |
---|
980 | endif |
---|
981 | |
---|
982 | IF (.NOT. child%var%interpolationshouldbemade) THEN |
---|
983 | C |
---|
984 | C |
---|
985 | C Calculation of the coefficients c1t and c2t for the temporary |
---|
986 | C interpolation |
---|
987 | if (pweight) then |
---|
988 | C |
---|
989 | c1t = weight |
---|
990 | C |
---|
991 | else |
---|
992 | C |
---|
993 | c1t = (REAL(AGRIF_Nbstepint()) + 1.) / Agrif_Rhot() |
---|
994 | C |
---|
995 | endif |
---|
996 | C |
---|
997 | c2t = 1. - c1t |
---|
998 | C |
---|
999 | C Time interpolation |
---|
1000 | C |
---|
1001 | kindex = 1 |
---|
1002 | C |
---|
1003 | do nb = 1,nbdim |
---|
1004 | C |
---|
1005 | do ndir = 1,2 |
---|
1006 | C |
---|
1007 | if (loctab_child(nb) /= (-ndir) |
---|
1008 | & .AND. loctab_child(nb) /= -3) then |
---|
1009 | |
---|
1010 | Call timeInterpolation |
---|
1011 | & (child,ptres2(:,:,ndir,nb),kindex,c1t,c2t,nbdim) |
---|
1012 | endif |
---|
1013 | C |
---|
1014 | enddo |
---|
1015 | C |
---|
1016 | enddo |
---|
1017 | C |
---|
1018 | |
---|
1019 | ENDIF |
---|
1020 | C |
---|
1021 | End Subroutine Agrif_Correctnd |
---|
1022 | C |
---|
1023 | C |
---|
1024 | C ************************************************************************** |
---|
1025 | CCC Subroutine saveAfterInterp |
---|
1026 | C ************************************************************************** |
---|
1027 | C |
---|
1028 | Subroutine saveAfterInterp(child,bounds,kindex,newsize,nbdim) |
---|
1029 | C |
---|
1030 | CCC Descritpion: |
---|
1031 | CCC Subroutine used to save the values of the grid variable on the fine grid |
---|
1032 | CCC after the space interpolation. |
---|
1033 | C |
---|
1034 | C Declarations: |
---|
1035 | C |
---|
1036 | |
---|
1037 | C |
---|
1038 | C argument |
---|
1039 | TYPE (AGRIF_PVariable) :: child ! The fine grid variable |
---|
1040 | INTEGER :: kindex ! Index indicating where this safeguard |
---|
1041 | ! is done on the fine grid |
---|
1042 | INTEGER :: nbdim, newsize |
---|
1043 | INTEGER,DIMENSION(nbdim,2) :: bounds |
---|
1044 | C |
---|
1045 | C Local scalars |
---|
1046 | INTEGER :: ir,jr,kr,lr,mr,nr |
---|
1047 | C |
---|
1048 | C |
---|
1049 | C Allocation of the array oldvalues2d |
---|
1050 | |
---|
1051 | C |
---|
1052 | if (newsize .LE. 0) return |
---|
1053 | C |
---|
1054 | Call Agrif_Checksize |
---|
1055 | & (child,kindex+newsize) |
---|
1056 | |
---|
1057 | if (child % var % interpIndex |
---|
1058 | & /= Agrif_Curgrid % parent % ngridstep ) then |
---|
1059 | child%var%oldvalues2d(1,kindex:kindex+newsize-1)= |
---|
1060 | & child%var%oldvalues2d(2,kindex:kindex+newsize-1) |
---|
1061 | endif |
---|
1062 | |
---|
1063 | SELECT CASE (nbdim) |
---|
1064 | CASE (1) |
---|
1065 | |
---|
1066 | !CDIR ALTCODE |
---|
1067 | do ir=bounds(1,1),bounds(1,2) |
---|
1068 | child%var%oldvalues2d(2,kindex) = |
---|
1069 | & child%var%parray1(ir) |
---|
1070 | kindex = kindex + 1 |
---|
1071 | enddo |
---|
1072 | C |
---|
1073 | CASE (2) |
---|
1074 | |
---|
1075 | do jr=bounds(2,1),bounds(2,2) |
---|
1076 | !CDIR ALTCODE |
---|
1077 | do ir=bounds(1,1),bounds(1,2) |
---|
1078 | child%var%oldvalues2d(2,kindex) = |
---|
1079 | & child%var%parray2(ir,jr) |
---|
1080 | kindex = kindex + 1 |
---|
1081 | enddo |
---|
1082 | enddo |
---|
1083 | C |
---|
1084 | CASE (3) |
---|
1085 | do kr=bounds(3,1),bounds(3,2) |
---|
1086 | do jr=bounds(2,1),bounds(2,2) |
---|
1087 | !CDIR ALTCODE |
---|
1088 | do ir=bounds(1,1),bounds(1,2) |
---|
1089 | child%var%oldvalues2d(2,kindex) = |
---|
1090 | & child%var%parray3(ir,jr,kr) |
---|
1091 | kindex = kindex + 1 |
---|
1092 | enddo |
---|
1093 | enddo |
---|
1094 | enddo |
---|
1095 | C |
---|
1096 | CASE (4) |
---|
1097 | do lr=bounds(4,1),bounds(4,2) |
---|
1098 | do kr=bounds(3,1),bounds(3,2) |
---|
1099 | do jr=bounds(2,1),bounds(2,2) |
---|
1100 | !CDIR ALTCODE |
---|
1101 | do ir=bounds(1,1),bounds(1,2) |
---|
1102 | child%var%oldvalues2d(2,kindex) = |
---|
1103 | & child%var%parray4(ir,jr,kr,lr) |
---|
1104 | kindex = kindex + 1 |
---|
1105 | enddo |
---|
1106 | enddo |
---|
1107 | enddo |
---|
1108 | enddo |
---|
1109 | C |
---|
1110 | CASE (5) |
---|
1111 | do mr=bounds(5,1),bounds(5,2) |
---|
1112 | do lr=bounds(4,1),bounds(4,2) |
---|
1113 | do kr=bounds(3,1),bounds(3,2) |
---|
1114 | do jr=bounds(2,1),bounds(2,2) |
---|
1115 | !CDIR ALTCODE |
---|
1116 | do ir=bounds(1,1),bounds(1,2) |
---|
1117 | child%var%oldvalues2d(2,kindex) = |
---|
1118 | & child%var%parray5(ir,jr,kr,lr,mr) |
---|
1119 | kindex = kindex + 1 |
---|
1120 | enddo |
---|
1121 | enddo |
---|
1122 | enddo |
---|
1123 | enddo |
---|
1124 | enddo |
---|
1125 | C |
---|
1126 | CASE (6) |
---|
1127 | do nr=bounds(6,1),bounds(6,2) |
---|
1128 | do mr=bounds(5,1),bounds(5,2) |
---|
1129 | do lr=bounds(4,1),bounds(4,2) |
---|
1130 | do kr=bounds(3,1),bounds(3,2) |
---|
1131 | do jr=bounds(2,1),bounds(2,2) |
---|
1132 | !CDIR ALTCODE |
---|
1133 | do ir=bounds(1,1),bounds(1,2) |
---|
1134 | child%var%oldvalues2d(2,kindex) = |
---|
1135 | & child%var%parray6(ir,jr,kr,lr,mr,nr) |
---|
1136 | kindex = kindex + 1 |
---|
1137 | enddo |
---|
1138 | enddo |
---|
1139 | enddo |
---|
1140 | enddo |
---|
1141 | enddo |
---|
1142 | enddo |
---|
1143 | END SELECT |
---|
1144 | C |
---|
1145 | C |
---|
1146 | End subroutine saveAfterInterp |
---|
1147 | C |
---|
1148 | C |
---|
1149 | C |
---|
1150 | C ************************************************************************** |
---|
1151 | CCC Subroutine timeInterpolation |
---|
1152 | C ************************************************************************** |
---|
1153 | C |
---|
1154 | Subroutine timeInterpolation(child,bounds,kindex,c1t,c2t,nbdim) |
---|
1155 | C |
---|
1156 | CCC Descritpion: |
---|
1157 | CCC Subroutine for a linear time interpolation on the child grid. |
---|
1158 | C |
---|
1159 | C Declarations: |
---|
1160 | C |
---|
1161 | |
---|
1162 | C |
---|
1163 | C argument |
---|
1164 | TYPE (AGRIF_PVariable) :: child ! The fine grid variable |
---|
1165 | INTEGER :: nbdim |
---|
1166 | INTEGER,DIMENSION(nbdim,2) :: bounds |
---|
1167 | INTEGER :: kindex ! Index indicating the values of the fine |
---|
1168 | ! grid got before and after the space |
---|
1169 | ! interpolation and used for the time |
---|
1170 | ! interpolation |
---|
1171 | REAL :: c1t,c2t! coefficients for the time interpolation |
---|
1172 | ! (c2t=1-c1t) |
---|
1173 | C |
---|
1174 | C Local aruments |
---|
1175 | INTEGER :: i |
---|
1176 | C Local scalars |
---|
1177 | INTEGER :: ir,jr,kr,lr,mr,nr |
---|
1178 | C |
---|
1179 | C |
---|
1180 | |
---|
1181 | SELECT CASE (nbdim) |
---|
1182 | CASE (1) |
---|
1183 | |
---|
1184 | !CDIR ALTCODE |
---|
1185 | do ir=bounds(1,1),bounds(1,2) |
---|
1186 | child%var%parray1(ir) = |
---|
1187 | & c2t*child % var % oldvalues2d(1,kindex) |
---|
1188 | & + c1t*child % var % oldvalues2d(2,kindex) |
---|
1189 | kindex = kindex + 1 |
---|
1190 | enddo |
---|
1191 | C |
---|
1192 | CASE (2) |
---|
1193 | |
---|
1194 | do jr=bounds(2,1),bounds(2,2) |
---|
1195 | !CDIR ALTCODE |
---|
1196 | do ir=bounds(1,1),bounds(1,2) |
---|
1197 | child%var%parray2(ir,jr) = |
---|
1198 | & c2t*child % var % oldvalues2d(1,kindex) |
---|
1199 | & + c1t*child % var % oldvalues2d(2,kindex) |
---|
1200 | kindex = kindex + 1 |
---|
1201 | enddo |
---|
1202 | enddo |
---|
1203 | C |
---|
1204 | CASE (3) |
---|
1205 | do kr=bounds(3,1),bounds(3,2) |
---|
1206 | do jr=bounds(2,1),bounds(2,2) |
---|
1207 | !CDIR ALTCODE |
---|
1208 | do ir=bounds(1,1),bounds(1,2) |
---|
1209 | child%var%parray3(ir,jr,kr) = |
---|
1210 | & c2t*child % var % oldvalues2d(1,kindex) |
---|
1211 | & + c1t*child % var % oldvalues2d(2,kindex) |
---|
1212 | kindex = kindex + 1 |
---|
1213 | enddo |
---|
1214 | enddo |
---|
1215 | enddo |
---|
1216 | C |
---|
1217 | CASE (4) |
---|
1218 | do lr=bounds(4,1),bounds(4,2) |
---|
1219 | do kr=bounds(3,1),bounds(3,2) |
---|
1220 | do jr=bounds(2,1),bounds(2,2) |
---|
1221 | !CDIR ALTCODE |
---|
1222 | do ir=bounds(1,1),bounds(1,2) |
---|
1223 | child%var%parray4(ir,jr,kr,lr) = |
---|
1224 | & c2t*child % var % oldvalues2d(1,kindex) |
---|
1225 | & + c1t*child % var % oldvalues2d(2,kindex) |
---|
1226 | kindex = kindex + 1 |
---|
1227 | enddo |
---|
1228 | enddo |
---|
1229 | enddo |
---|
1230 | enddo |
---|
1231 | C |
---|
1232 | CASE (5) |
---|
1233 | do mr=bounds(5,1),bounds(5,2) |
---|
1234 | do lr=bounds(4,1),bounds(4,2) |
---|
1235 | do kr=bounds(3,1),bounds(3,2) |
---|
1236 | do jr=bounds(2,1),bounds(2,2) |
---|
1237 | !CDIR ALTCODE |
---|
1238 | do ir=bounds(1,1),bounds(1,2) |
---|
1239 | child%var%parray5(ir,jr,kr,lr,mr) = |
---|
1240 | & c2t*child % var % oldvalues2d(1,kindex) |
---|
1241 | & + c1t*child % var % oldvalues2d(2,kindex) |
---|
1242 | kindex = kindex + 1 |
---|
1243 | enddo |
---|
1244 | enddo |
---|
1245 | enddo |
---|
1246 | enddo |
---|
1247 | enddo |
---|
1248 | C |
---|
1249 | CASE (6) |
---|
1250 | do nr=bounds(6,1),bounds(6,2) |
---|
1251 | do mr=bounds(5,1),bounds(5,2) |
---|
1252 | do lr=bounds(4,1),bounds(4,2) |
---|
1253 | do kr=bounds(3,1),bounds(3,2) |
---|
1254 | do jr=bounds(2,1),bounds(2,2) |
---|
1255 | !CDIR ALTCODE |
---|
1256 | do ir=bounds(1,1),bounds(1,2) |
---|
1257 | child%var%parray6(ir,jr,kr,lr,mr,nr) = |
---|
1258 | & c2t*child % var % oldvalues2d(1,kindex) |
---|
1259 | & + c1t*child % var % oldvalues2d(2,kindex) |
---|
1260 | kindex = kindex + 1 |
---|
1261 | enddo |
---|
1262 | enddo |
---|
1263 | enddo |
---|
1264 | enddo |
---|
1265 | enddo |
---|
1266 | enddo |
---|
1267 | END SELECT |
---|
1268 | |
---|
1269 | C |
---|
1270 | C |
---|
1271 | End subroutine timeInterpolation |
---|
1272 | C |
---|
1273 | C |
---|
1274 | C |
---|
1275 | C ************************************************************************** |
---|
1276 | CCC Subroutine Agrif_Checksize |
---|
1277 | C ************************************************************************** |
---|
1278 | C |
---|
1279 | Subroutine Agrif_Checksize(child,newsize) |
---|
1280 | C |
---|
1281 | CCC Descritpion: |
---|
1282 | CCC Subroutine used in the saveAfterInterp procedure to allocate the |
---|
1283 | CCC oldvalues2d array. |
---|
1284 | C |
---|
1285 | C Declarations: |
---|
1286 | C |
---|
1287 | |
---|
1288 | C |
---|
1289 | C TYPE argument |
---|
1290 | TYPE (AGRIF_PVariable) :: child ! The fine grid variable |
---|
1291 | C |
---|
1292 | C Scalar arguments |
---|
1293 | INTEGER :: newsize ! Size of the domains where the boundary |
---|
1294 | ! conditions are calculated |
---|
1295 | C |
---|
1296 | C Local arrays |
---|
1297 | REAL, DIMENSION(:,:), Allocatable :: tempoldvalues ! Temporary array |
---|
1298 | C |
---|
1299 | C |
---|
1300 | if (.NOT. associated(child % var % oldvalues2d)) then |
---|
1301 | C |
---|
1302 | allocate(child % var % oldvalues2d(2,newsize)) |
---|
1303 | C |
---|
1304 | child % var % oldvalues2d=0. |
---|
1305 | C |
---|
1306 | else |
---|
1307 | C |
---|
1308 | if (SIZE(child % var % oldvalues2d,2) < newsize) then |
---|
1309 | C |
---|
1310 | allocate(tempoldvalues(2,SIZE(child % var % |
---|
1311 | & oldvalues2d,2))) |
---|
1312 | C |
---|
1313 | tempoldvalues = child % var % oldvalues2d |
---|
1314 | C |
---|
1315 | deallocate(child % var % oldvalues2d) |
---|
1316 | C |
---|
1317 | allocate(child % var % oldvalues2d(2,newsize)) |
---|
1318 | C |
---|
1319 | child%var%oldvalues2d=0. |
---|
1320 | C |
---|
1321 | child % var % oldvalues2d(:,1:SIZE(tempoldvalues,2)) = |
---|
1322 | & tempoldvalues(:,:) |
---|
1323 | C |
---|
1324 | deallocate(tempoldvalues) |
---|
1325 | C |
---|
1326 | endif |
---|
1327 | C |
---|
1328 | endif |
---|
1329 | C |
---|
1330 | C |
---|
1331 | End Subroutine Agrif_Checksize |
---|
1332 | C |
---|
1333 | C |
---|
1334 | C |
---|
1335 | C |
---|
1336 | End Module AGRIF_boundary |
---|
1337 | |
---|