New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
modmask.F90 in branches/UKMO/dev_r5518_GO6_package/NEMOGCM/EXTERNAL/AGRIF/AGRIF_FILES – NEMO

source: branches/UKMO/dev_r5518_GO6_package/NEMOGCM/EXTERNAL/AGRIF/AGRIF_FILES/modmask.F90 @ 7993

Last change on this file since 7993 was 7993, checked in by frrh, 7 years ago

Merge in missing revisions 6428:2477 inclusive and 6482 from nemo_v3_6_STABLE
branch. In ptic, this includes the fix for restartability of runoff fields in coupled
models. Evolution of coupled models will therefor be affected.

These changes donot affect evolution of the current stand-alone NEMO-CICE GO6
standard configuration.

Work and testing documented in Met Office GMED ticket 320.

File size: 20.1 KB
Line 
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_Mask.
25!>
26!> Module for masks.
27!
28module Agrif_Mask
29!
30    use Agrif_Types
31!
32    implicit none
33!
34contains
35!
36!===================================================================================================
37!  subroutine Agrif_CheckMasknD
38!
39!> Called in the procedure #Agrif_InterpnD to recalculate the value of the parent grid variable
40!! when this one is equal to Agrif_SpecialValue.
41!---------------------------------------------------------------------------------------------------
42subroutine Agrif_CheckMasknD ( tempP, parent, pbtab, petab, ppbtab, ppetab, noraftab, nbdim )
43!---------------------------------------------------------------------------------------------------
44    type(Agrif_Variable), pointer   :: tempP  !< Part of the parent grid used for the interpolation of the child grid
45    type(Agrif_Variable), pointer   :: parent !< The parent grid
46    integer, dimension(nbdim)       :: pbtab  !< limits of the parent grid used
47    integer, dimension(nbdim)       :: petab  !< interpolation of the child grid
48    integer, dimension(nbdim)       :: ppbtab, ppetab
49    logical, dimension(nbdim)       :: noraftab
50    integer                         :: nbdim
51!
52    integer :: i0,j0,k0,l0,m0,n0
53!
54    select case (nbdim)
55    case (1)
56        do i0 = pbtab(1),petab(1)
57            if (tempP%array1(i0) == Agrif_SpecialValue) then
58                call CalculNewValTempP((/i0/),tempP,parent,ppbtab,ppetab,noraftab,nbdim)
59            endif
60        enddo
61    case (2)
62        do j0 = pbtab(2),petab(2)
63        do i0 = pbtab(1),petab(1)
64            if (tempP%array2(i0,j0) == Agrif_SpecialValue) then
65                call CalculNewValTempP((/i0,j0/),tempP,parent,ppbtab,ppetab, noraftab,nbdim)
66            endif
67        enddo
68        enddo
69    case (3)
70        do k0 = pbtab(3),petab(3)
71        do j0 = pbtab(2),petab(2)
72        do i0 = pbtab(1),petab(1)
73            if (tempP%array3(i0,j0,k0) == Agrif_SpecialValue) then
74!------CDIR NEXPAND
75                call CalculNewValTempP3D((/i0,j0,k0/), &
76                    tempP%array3(ppbtab(1),ppbtab(2),ppbtab(3)),  &
77                    parent%array3(ppbtab(1),ppbtab(2),ppbtab(3)), &
78                    ppbtab,ppetab,noraftab,MaxSearch,Agrif_SpecialValue)
79
80!            Call CalculNewValTempP((/i0,j0,k0/),
81!     &                             tempP,parent,
82!     &                             ppbtab,ppetab,
83!     &                             noraftab,nbdim)
84
85            endif
86        enddo
87        enddo
88        enddo
89    case (4)
90        do l0 = pbtab(4),petab(4)
91        do k0 = pbtab(3),petab(3)
92        do j0 = pbtab(2),petab(2)
93        do i0 = pbtab(1),petab(1)
94            if (tempP%array4(i0,j0,k0,l0) == Agrif_SpecialValue) then
95                call CalculNewValTempP4D((/i0,j0,k0,l0/), &
96                    tempP%array4(ppbtab(1),ppbtab(2),ppbtab(3),ppbtab(4)),  &
97                    parent%array4(ppbtab(1),ppbtab(2),ppbtab(3),ppbtab(4)), &
98                    ppbtab,ppetab,noraftab,MaxSearch,Agrif_SpecialValue)
99            endif
100        enddo
101        enddo
102        enddo
103        enddo
104    case (5)
105        do m0 = pbtab(5),petab(5)
106        do l0 = pbtab(4),petab(4)
107        do k0 = pbtab(3),petab(3)
108        do j0 = pbtab(2),petab(2)
109        do i0 = pbtab(1),petab(1)
110            if (tempP%array5(i0,j0,k0,l0,m0) == Agrif_SpecialValue) then
111                call CalculNewValTempP((/i0,j0,k0,l0,m0/), &
112                    tempP,parent,ppbtab,ppetab,noraftab,nbdim)
113            endif
114        enddo
115        enddo
116        enddo
117        enddo
118        enddo
119    case (6)
120        do n0 = pbtab(6),petab(6)
121        do m0 = pbtab(5),petab(5)
122        do l0 = pbtab(4),petab(4)
123        do k0 = pbtab(3),petab(3)
124        do j0 = pbtab(2),petab(2)
125        do i0 = pbtab(1),petab(1)
126            if (tempP%array6(i0,j0,k0,l0,m0,n0) == Agrif_SpecialValue) then
127                call CalculNewValTempP((/i0,j0,k0,l0,m0,n0/), &
128                    tempP,parent,ppbtab,ppetab,noraftab,nbdim)
129            endif
130        enddo
131        enddo
132        enddo
133        enddo
134        enddo
135        enddo
136    end select
137!---------------------------------------------------------------------------------------------------
138end subroutine Agrif_CheckMasknD
139!===================================================================================================
140!
141!===================================================================================================
142!  subroutine CalculNewValTempP
143!
144!> Called in the procedure #Agrif_InterpnD to recalculate the value of the parent grid variable
145!! when this one is equal to Agrif_SpecialValue.
146!---------------------------------------------------------------------------------------------------
147subroutine CalculNewValTempP ( indic, tempP, parent, ppbtab, ppetab, noraftab, nbdim )
148!---------------------------------------------------------------------------------------------------
149    integer, dimension(nbdim)       :: indic
150    type(Agrif_Variable), pointer   :: tempP  !< Part of the parent grid used for the interpolation of the child grid
151    type(Agrif_Variable), pointer   :: parent !< The parent grid
152    integer, dimension(nbdim)       :: ppbtab, ppetab
153    logical, dimension(nbdim)       :: noraftab
154    integer                         :: nbdim
155!
156    integer                     :: i,ii,iii,jj,kk,ll,mm,nn
157    integer, dimension(nbdim)   :: imin,imax,idecal
158    integer                     :: nbvals
159    real                        :: res
160    real                        :: valparent
161    integer                     :: ValMax
162    logical                     :: firsttest
163!
164    ValMax = 1
165    do iii = 1,nbdim
166        if (.NOT.noraftab(iii)) then
167            ValMax = max(ValMax,ppetab(iii)-indic(iii))
168            ValMax = max(ValMax,indic(iii)-ppbtab(iii))
169        endif
170    enddo
171!
172    Valmax = min(Valmax,MaxSearch)
173!
174!CDIR NOVECTOR
175    imin = indic
176!CDIR NOVECTOR
177    imax = indic
178!
179    i = 1
180    firsttest = .TRUE.
181!
182    do while (i <= ValMax)
183!
184        if ( (i == 1).AND.(firsttest) ) i = Valmax
185!
186        do iii = 1,nbdim
187            if (.NOT.noraftab(iii)) then
188                imin(iii) = max(indic(iii) - i,ppbtab(iii))
189                imax(iii) = min(indic(iii) + i,ppetab(iii))
190                if (firsttest) then
191                    if (indic(iii) > ppbtab(iii)) then
192!CDIR NOVECTOR
193                        idecal = indic
194                        idecal(iii) = idecal(iii)-1
195                        SELECT CASE(nbdim)
196                        CASE (1)
197                            if (tempP%array1(idecal(1) &
198                                    ) == Agrif_SpecialValue) imin(iii) = imax(iii)
199                        CASE (2)
200                            if (tempP%array2(idecal(1), idecal(2) &
201                                    ) == Agrif_SpecialValue) imin(iii) = imax(iii)
202                        CASE (3)
203                            if (tempP%array3(idecal(1), &
204                                             idecal(2), idecal(3) &
205                                    ) == Agrif_SpecialValue) imin(iii) = imax(iii)
206                        CASE (4)
207                            if (tempP%array4(idecal(1), idecal(2), &
208                                             idecal(3), idecal(4)  &
209                                    ) == Agrif_SpecialValue) imin(iii) = imax(iii)
210                        CASE (5)
211                            if (tempP%array5(idecal(1), idecal(2), &
212                                             idecal(3), idecal(4), &
213                                             idecal(5)             &
214                                    ) == Agrif_SpecialValue) imin(iii) = imax(iii)
215                        CASE (6)
216                            if (tempP%array6(idecal(1), idecal(2), &
217                                             idecal(3), idecal(4), &
218                                             idecal(5), idecal(6)  &
219                                    ) == Agrif_SpecialValue) imin(iii) = imax(iii)
220                        END SELECT
221                    endif
222                endif
223            endif
224        enddo
225!
226        Res = 0.
227        Nbvals = 0
228!
229        SELECT CASE(nbdim)
230        CASE (1)
231!CDIR ALTCODE
232!CDIR SHORTLOOP
233            do ii = imin(1),imax(1)
234                ValParent = parent%array1(ii)
235                if ( ValParent /= Agrif_SpecialValue) then
236                    Res = Res + ValParent
237                    Nbvals = Nbvals + 1
238                endif
239            enddo
240!
241        CASE (2)
242            do jj = imin(2),imax(2)
243!CDIR ALTCODE
244!CDIR SHORTLOOP
245            do ii = imin(1),imax(1)
246                ValParent = parent%array2(ii,jj)
247                if ( ValParent /= Agrif_SpecialValue) then
248                    Res = Res + ValParent
249                    Nbvals = Nbvals + 1
250                endif
251            enddo
252            enddo
253
254        CASE (3)
255            do kk = imin(3),imax(3)
256            do jj = imin(2),imax(2)
257!CDIR ALTCODE
258!CDIR SHORTLOOP
259            do ii = imin(1),imax(1)
260                ValParent = parent%array3(ii,jj,kk)
261                if ( ValParent /= Agrif_SpecialValue) then
262                    Res = Res + ValParent
263                    Nbvals = Nbvals + 1
264                endif
265            enddo
266            enddo
267            enddo
268
269        CASE (4)
270            do ll = imin(4),imax(4)
271            do kk = imin(3),imax(3)
272            do jj = imin(2),imax(2)
273!CDIR ALTCODE
274!CDIR SHORTLOOP
275            do ii = imin(1),imax(1)
276                ValParent = parent%array4(ii,jj,kk,ll)
277                if ( ValParent /= Agrif_SpecialValue) then
278                    Res = Res + ValParent
279                    Nbvals = Nbvals + 1
280                endif
281            enddo
282            enddo
283            enddo
284            enddo
285
286        CASE (5)
287            do mm = imin(5),imax(5)
288            do ll = imin(4),imax(4)
289            do kk = imin(3),imax(3)
290            do jj = imin(2),imax(2)
291!CDIR ALTCODE
292!CDIR SHORTLOOP
293            do ii = imin(1),imax(1)
294                ValParent = parent%array5(ii,jj,kk,ll,mm)
295                if ( ValParent /= Agrif_SpecialValue) then
296                    Res = Res + ValParent
297                    Nbvals = Nbvals + 1
298                endif
299            enddo
300            enddo
301            enddo
302            enddo
303            enddo
304
305        CASE (6)
306            do nn = imin(6),imax(6)
307            do mm = imin(5),imax(5)
308            do ll = imin(4),imax(4)
309            do kk = imin(3),imax(3)
310            do jj = imin(2),imax(2)
311!CDIR ALTCODE
312!CDIR SHORTLOOP
313            do ii = imin(1),imax(1)
314                ValParent = parent%array6(ii,jj,kk,ll,mm,nn)
315                if ( ValParent /= Agrif_SpecialValue) then
316                    Res = Res + ValParent
317                    Nbvals = Nbvals + 1
318                endif
319            enddo
320            enddo
321            enddo
322            enddo
323            enddo
324            enddo
325
326        END SELECT
327!
328        if (Nbvals > 0) then
329            if (firsttest) then
330                firsttest = .FALSE.
331                i=1
332                cycle
333            endif
334            SELECT CASE(nbdim)
335            CASE (1)
336                tempP%array1(indic(1))           = Res/Nbvals
337            CASE (2)
338                tempP%array2(indic(1), indic(2)) = Res/Nbvals
339            CASE (3)
340                tempP%array3(indic(1), indic(2), &
341                             indic(3))           = Res/Nbvals
342            CASE (4)
343                tempP%array4(indic(1), indic(2), &
344                             indic(3), indic(4)) = Res/Nbvals
345            CASE (5)
346                tempP%array5(indic(1), indic(2), &
347                             indic(3), indic(4), &
348                             indic(5))           = Res/Nbvals
349            CASE (6)
350                tempP%array6(indic(1), indic(2), &
351                             indic(3), indic(4), &
352                             indic(5), indic(6)) = Res/Nbvals
353            END SELECT
354            exit
355        else
356            if (firsttest) exit
357            i = i + 1
358        endif
359!
360    enddo
361!---------------------------------------------------------------------------------------------------
362end subroutine CalculNewValTempP
363!===================================================================================================
364!
365!===================================================================================================
366!  subroutine CalculNewValTempP3D
367!
368!> Called in the procedure #Agrif_InterpnD to recalculate the value of the parent grid variable
369!! when this one is equal to Agrif_SpecialValue.
370!---------------------------------------------------------------------------------------------------
371subroutine CalculNewValTempP3D ( indic, tempP, parent, ppbtab, ppetab, noraftab, &
372                                 MaxSearch, Agrif_SpecialValue )
373!---------------------------------------------------------------------------------------------------
374    integer, parameter          :: nbdim = 3
375    integer, dimension(nbdim)   :: indic
376    integer, dimension(nbdim)   :: ppbtab, ppetab
377    logical, dimension(nbdim)   :: noraftab
378    integer                     :: MaxSearch
379    real                        :: Agrif_SpecialValue
380    real, dimension(ppbtab(1):ppetab(1), &
381                    ppbtab(2):ppetab(2), &
382                    ppbtab(3):ppetab(3)) &
383                                :: tempP, parent  !< Part of the parent grid used for
384                                                  !< the interpolation of the child grid
385!
386    integer                     :: i,ii,iii,jj,kk
387    integer, dimension(nbdim)   :: imin,imax,idecal
388    integer                     :: Nbvals
389    real                        :: Res
390    real                        :: ValParent
391    integer                     :: ValMax
392    logical                     :: Existunmasked
393!
394    ValMax = 1
395!CDIR NOVECTOR
396    do iii = 1,nbdim
397        if (.NOT.noraftab(iii)) then
398            ValMax = max(ValMax,ppetab(iii)-indic(iii))
399            ValMax = max(ValMax,indic(iii)-ppbtab(iii))
400        endif
401    enddo
402!
403    Valmax = min(Valmax,MaxSearch)
404!
405!CDIR NOVECTOR
406    imin = indic
407!CDIR NOVECTOR
408    imax = indic
409
410!CDIR NOVECTOR
411    idecal = indic
412    i = Valmax
413!
414    do iii = 1,nbdim
415        if (.NOT.noraftab(iii)) then
416            imin(iii) = max(indic(iii) - i,ppbtab(iii))
417            imax(iii) = min(indic(iii) + i,ppetab(iii))
418
419            if (indic(iii) > ppbtab(iii)) then
420                idecal(iii) = idecal(iii)-1
421                if (tempP(idecal(1),idecal(2),idecal(3)) == Agrif_SpecialValue) then
422                    imin(iii) = imax(iii)
423                endif
424                idecal(iii) = idecal(iii)+1
425            endif
426        endif
427    enddo
428!
429    Existunmasked = .FALSE.
430!
431    do kk = imin(3),imax(3)
432    do jj = imin(2),imax(2)
433!CDIR NOVECTOR
434    do ii = imin(1),imax(1)
435        if ( parent(ii,jj,kk) /= Agrif_SpecialValue) then
436            Existunmasked = .TRUE.
437            exit
438        endif
439    enddo
440    enddo
441    enddo
442!
443    if (.Not.Existunmasked) return
444!
445    i = 1
446!
447    do while(i <= ValMax)
448!
449        do iii = 1 , nbdim
450            if (.NOT.noraftab(iii)) then
451                imin(iii) = max(indic(iii) - i,ppbtab(iii))
452                imax(iii) = min(indic(iii) + i,ppetab(iii))
453            endif
454        enddo
455!
456        Res = 0.
457        Nbvals = 0
458!
459        do kk = imin(3),imax(3)
460        do jj = imin(2),imax(2)
461!CDIR NOVECTOR
462        do ii = imin(1),imax(1)
463            ValParent = parent(ii,jj,kk)
464            if ( ValParent /= Agrif_SpecialValue) then
465                Res = Res + ValParent
466                Nbvals = Nbvals + 1
467            endif
468        enddo
469        enddo
470        enddo
471!
472        if (Nbvals > 0) then
473            tempP(indic(1),indic(2),indic(3)) = Res/Nbvals
474            exit
475        else
476            i = i + 1
477        endif
478    enddo
479!---------------------------------------------------------------------------------------------------
480end subroutine CalculNewValTempP3D
481!===================================================================================================
482!
483!===================================================================================================
484!  subroutine CalculNewValTempP4D
485!
486!> Called in the procedure #Agrif_InterpnD to recalculate the value of the parent grid variable
487!! when this one is equal to Agrif_SpecialValue.
488!---------------------------------------------------------------------------------------------------
489subroutine CalculNewValTempP4D ( indic, tempP, parent, ppbtab, ppetab, noraftab, &
490                                 MaxSearch, Agrif_SpecialValue )
491!---------------------------------------------------------------------------------------------------
492    integer, parameter          :: nbdim = 4
493    integer, dimension(nbdim)   :: indic
494    integer, dimension(nbdim)   :: ppbtab, ppetab
495    logical, dimension(nbdim)   :: noraftab
496    integer                     :: MaxSearch
497    real                        :: Agrif_SpecialValue
498    real, dimension(ppbtab(1):ppetab(1), &
499                    ppbtab(2):ppetab(2), &
500                    ppbtab(3):ppetab(3), &
501                    ppbtab(4):ppetab(4)) &
502                                :: tempP, parent  !< Part of the parent grid used for
503                                                  !< the interpolation of the child grid
504!
505    integer                     :: i,ii,iii,jj,kk,ll
506    integer, dimension(nbdim)   :: imin,imax,idecal
507    integer                     :: Nbvals
508    real                        :: Res
509    real                        :: ValParent
510    integer                     :: ValMax
511!
512    logical                     :: firsttest
513!
514    ValMax = 1
515    do iii = 1,nbdim
516        if (.NOT.noraftab(iii)) then
517            ValMax = max(ValMax,ppetab(iii)-indic(iii))
518            ValMax = max(ValMax,indic(iii)-ppbtab(iii))
519        endif
520    enddo
521!
522    Valmax = min(Valmax,MaxSearch)
523!
524    imin = indic
525    imax = indic
526!
527    i = 1
528    firsttest = .TRUE.
529    idecal = indic
530!
531    do while (i <= ValMax)
532!
533        if ((i == 1).AND.(firsttest)) i = Valmax
534
535        do iii = 1,nbdim
536            if (.NOT.noraftab(iii)) then
537                imin(iii) = max(indic(iii) - i,ppbtab(iii))
538                imax(iii) = min(indic(iii) + i,ppetab(iii))
539                if (firsttest) then
540                    if (indic(iii) > ppbtab(iii)) then
541                        idecal(iii) = idecal(iii)-1
542                        if (tempP(idecal(1),idecal(2),idecal(3),idecal(4)) == Agrif_SpecialValue) then
543                            imin(iii) = imax(iii)
544                        endif
545                        idecal(iii) = idecal(iii)+1
546                    endif
547                endif
548            endif
549        enddo
550!
551        Res = 0.
552        Nbvals = 0
553!
554        do ll = imin(4),imax(4)
555        do kk = imin(3),imax(3)
556        do jj = imin(2),imax(2)
557        do ii = imin(1),imax(1)
558            ValParent = parent(ii,jj,kk,ll)
559            if ( ValParent /= Agrif_SpecialValue) then
560                Res = Res + ValParent
561                Nbvals = Nbvals + 1
562            endif
563        enddo
564        enddo
565        enddo
566        enddo
567!
568        if (Nbvals > 0) then
569            if (firsttest) then
570                firsttest = .FALSE.
571                i=1
572                cycle
573            endif
574
575            tempP(indic(1),indic(2),indic(3),indic(4)) = Res/Nbvals
576            exit
577        else
578            if (firsttest) exit
579            i = i + 1
580        endif
581    enddo
582!---------------------------------------------------------------------------------------------------
583end subroutine CalculNewValTempP4D
584!===================================================================================================
585!
586end module Agrif_Mask
Note: See TracBrowser for help on using the repository browser.