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.F in branches/DEV_r2006_merge_TRA_TRC/AGRIF/AGRIF_FILES – NEMO

source: branches/DEV_r2006_merge_TRA_TRC/AGRIF/AGRIF_FILES/modmask.F @ 3097

Last change on this file since 3097 was 1200, checked in by rblod, 16 years ago

Adapt Agrif to the new SBC and correct several bugs for agrif (restart writing and reading), see ticket #133
Note : this fix does not work yet on NEC computerq (sxf90/360)

  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 21.8 KB
Line 
1!
2! $Id$
3!
4C     AGRIF (Adaptive Grid Refinement In Fortran)
5C
6C     Copyright (C) 2003 Laurent Debreu (Laurent.Debreu@imag.fr)
7C                        Christophe Vouland (Christophe.Vouland@imag.fr)   
8C
9C     This program is free software; you can redistribute it and/or modify
10C     it under the terms of the GNU General Public License as published by
11C     the Free Software Foundation; either version 2 of the License, or
12C     (at your option) any later version.
13C
14C     This program is distributed in the hope that it will be useful,
15C     but WITHOUT ANY WARRANTY; without even the implied warranty of
16C     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
17C     GNU General Public License for more details.
18C
19C     You should have received a copy of the GNU General Public License
20C     along with this program; if not, write to the Free Software
21C     Foundation, Inc., 59 Temple Place -  Suite 330, Boston, MA 02111-1307, USA.
22C
23C
24C
25CCC   Module Agrif_Mask
26C
27      Module Agrif_Mask
28C
29CCC   Description:
30CCC   Module for masks
31C
32C     Modules used: 
33C
34      Use Agrif_Types       
35C
36      IMPLICIT NONE
37C
38      CONTAINS
39C     Define procedures contained in this module
40C
41C     **************************************************************************
42C     Subroutine Agrif_CheckMasknD
43C     **************************************************************************
44C       
45      Subroutine Agrif_CheckMasknD(tempP,parent,pbtab,petab,ppbtab,
46     &               ppetab,noraftab,nbdim)
47C
48CCC   Description:
49CCC   Subroutine called in the procedure Agrif_InterpnD to recalculate the value
50CCC   of the parent grid variable when this one is equal to Agrif_SpecialValue. 
51C
52C     Declarations:
53C
54       
55C
56C     Arrays arguments     
57      INTEGER :: nbdim
58      INTEGER,DIMENSION(nbdim) :: pbtab  ! Limits of the parent grid used 
59      INTEGER,DIMENSION(nbdim) :: petab  ! interpolation of the child grid
60      LOGICAL,DIMENSION(nbdim) :: noraftab
61      INTEGER,DIMENSION(nbdim) :: ppbtab,ppetab
62C
63C     Pointer argument
64      TYPE(AGRIF_PVARIABLE) :: tempP  ! Part of the parent grid used for
65                                      ! the interpolation of the child grid                     
66C     Data TYPE argument                                   
67      TYPE(AGRIF_PVARIABLE) :: parent      ! The parent grid
68C
69C     Local scalar
70      INTEGER                   :: i0,j0,k0,l0,m0,n0
71C     
72C     Local arrays
73C
74C     
75      SELECT CASE (nbdim)
76      CASE (1)
77         do i0 = pbtab(1),petab(1)
78         IF (tempP%var%array1(i0)
79     &                        == Agrif_SpecialValue) Then
80            Call CalculNewValTempP((/i0/),
81     &                             tempP,parent,
82     &                             ppbtab,ppetab,
83     &                             noraftab,nbdim)
84         ENDIF
85         enddo
86      CASE (2)
87         do j0 = pbtab(2),petab(2)
88         do i0 = pbtab(1),petab(1)
89         IF (tempP%var%array2(i0,j0)
90     &                        == Agrif_SpecialValue) Then
91            Call CalculNewValTempP((/i0,j0/),
92     &                             tempP,parent,
93     &                             ppbtab,ppetab,
94     &                             noraftab,nbdim)
95         ENDIF
96         enddo 
97         enddo
98      CASE (3)
99         do k0 = pbtab(3),petab(3)
100         do j0 = pbtab(2),petab(2)
101         do i0 = pbtab(1),petab(1)
102         IF (tempP%var%array3(i0,j0,k0)
103     &                        == Agrif_SpecialValue) Then
104!------CDIR NEXPAND
105            Call CalculNewValTempP3D((/i0,j0,k0/),
106     &      tempP%var%array3(ppbtab(1),ppbtab(2),ppbtab(3)),
107     &      parent%var%array3(ppbtab(1),ppbtab(2),ppbtab(3)),
108     &                           ppbtab,ppetab,
109     &                           noraftab,MaxSearch,Agrif_SpecialValue)
110     
111c            Call CalculNewValTempP((/i0,j0,k0/),
112c     &                             tempP,parent,
113c     &                             ppbtab,ppetab,
114c     &                             noraftab,nbdim)
115         
116         ENDIF
117         enddo
118         enddo 
119         enddo
120      CASE (4)
121         do l0 = pbtab(4),petab(4)
122         do k0 = pbtab(3),petab(3)
123         do j0 = pbtab(2),petab(2)
124         do i0 = pbtab(1),petab(1)
125         IF (tempP%var%array4(i0,j0,k0,l0) 
126     &                        == Agrif_SpecialValue) Then
127            Call CalculNewValTempP4D((/i0,j0,k0,l0/),
128     &      tempP%var%array4(ppbtab(1),ppbtab(2),ppbtab(3),ppbtab(4)),
129     &      parent%var%array4(ppbtab(1),ppbtab(2),ppbtab(3),ppbtab(4)),
130     &                           ppbtab,ppetab,
131     &                           noraftab,MaxSearch,Agrif_SpecialValue)
132         ENDIF
133         enddo
134         enddo
135         enddo 
136         enddo
137      CASE (5)
138         do m0 = pbtab(5),petab(5)
139         do l0 = pbtab(4),petab(4)
140         do k0 = pbtab(3),petab(3)
141         do j0 = pbtab(2),petab(2)
142         do i0 = pbtab(1),petab(1)
143         IF (tempP%var%array5(i0,j0,k0,l0,m0) 
144     &                       == Agrif_SpecialValue) Then
145            Call CalculNewValTempP((/i0,j0,k0,l0,m0/),
146     &                             tempP,parent,
147     &                             ppbtab,ppetab,
148     &                             noraftab,nbdim)
149         ENDIF
150         enddo
151         enddo
152         enddo 
153         enddo
154         enddo
155      CASE (6)
156         do n0 = pbtab(6),petab(6)
157         do m0 = pbtab(5),petab(5)
158         do l0 = pbtab(4),petab(4)
159         do k0 = pbtab(3),petab(3)
160         do j0 = pbtab(2),petab(2)
161         do i0 = pbtab(1),petab(1)
162         IF (tempP%var%array6(i0,j0,k0,l0,m0,n0) 
163     &                       == Agrif_SpecialValue) Then
164            Call CalculNewValTempP((/i0,j0,k0,l0,m0,n0/),
165     &                             tempP,parent,
166     &                             ppbtab,ppetab,
167     &                             noraftab,nbdim)
168         ENDIF
169         enddo
170         enddo
171         enddo
172         enddo 
173         enddo
174         enddo
175      END SELECT
176C
177C     
178C     
179      End Subroutine Agrif_CheckMasknD
180C
181C
182C     **************************************************************************
183C     Subroutine CalculNewValTempP
184C     **************************************************************************
185C       
186      Subroutine CalculNewValTempP(indic,
187     &               tempP,parent,ppbtab,
188     &               ppetab,noraftab,nbdim)
189C
190CCC   Description:
191CCC   Subroutine called in the procedure Agrif_InterpnD to recalculate the value
192CCC   of the parent grid variable when this one is equal to Agrif_SpecialValue. 
193C
194C     Declarations:
195C
196       
197C
198C     Arrays arguments     
199      INTEGER :: nbdim
200      INTEGER,DIMENSION(nbdim) :: indic
201      LOGICAL,DIMENSION(nbdim) :: noraftab
202      INTEGER,DIMENSION(nbdim) :: ppbtab,ppetab
203C
204C     Pointer argument
205      TYPE(AGRIF_PVARIABLE) :: tempP  ! Part of the parent grid used for
206                                      ! the interpolation of the child grid                     
207C     Data TYPE argument                                   
208      TYPE(AGRIF_PVARIABLE) :: parent      ! The parent grid
209C
210C     Local scalar
211      INTEGER                  :: i,ii,iii,jj,kk,ll,mm,nn 
212      INTEGER,DIMENSION(nbdim) :: imin,imax,idecal
213      INTEGER                  :: Nbvals
214      REAL                     :: Res
215      REAL                     :: ValParent
216      INTEGER                  :: ValMax
217      LOGICAL                  :: firsttest
218C     
219C     Local arrays     
220C
221      ValMax = 1
222      do iii = 1 , nbdim
223         IF (.NOT.noraftab(iii)) THEN
224         ValMax = max(ValMax,ppetab(iii)-indic(iii))
225         ValMax = max(ValMax,indic(iii)-ppbtab(iii))
226         ENDIF
227      enddo
228C
229      Valmax = min(Valmax,MaxSearch)
230C
231!CDIR NOVECTOR
232      imin = indic
233!CDIR NOVECTOR
234      imax = indic
235C
236         i = 1
237         firsttest = .TRUE.
238C
239         do While(i <= ValMax)
240C
241         IF ((i == 1).AND.(firsttest)) i = Valmax
242
243            do iii = 1 , nbdim
244               if (.NOT.noraftab(iii)) then
245                  imin(iii) = max(indic(iii) - i,ppbtab(iii))
246                  imax(iii) = min(indic(iii) + i,ppetab(iii))
247                  if (firsttest) then               
248                  if (indic(iii).GT.ppbtab(iii)) then
249
250!CDIR NOVECTOR
251                     idecal = indic
252                     idecal(iii) = idecal(iii)-1
253                     SELECT CASE(nbdim)
254                     CASE (1)
255                        if (tempP%var%array1(idecal(1)
256     &                            ) == Agrif_SpecialValue) then
257                           imin(iii) = imax(iii)
258                        endif                     
259                     CASE (2)
260                        if (tempP%var%array2(idecal(1),
261     &            idecal(2)) == Agrif_SpecialValue) then
262                           imin(iii) = imax(iii)
263                        endif
264                     CASE (3)
265                        if (tempP%var%array3(idecal(1),
266     &            idecal(2),idecal(3)) 
267     &                               == Agrif_SpecialValue) then
268                           imin(iii) = imax(iii)
269                        endif 
270                     CASE (4)
271                        if (tempP%var%array4(idecal(1),
272     &            idecal(2),idecal(3),idecal(4)) 
273     &                               == Agrif_SpecialValue) then
274                           imin(iii) = imax(iii)
275                        endif     
276                     CASE (5)
277                        if (tempP%var%array5(idecal(1),
278     &            idecal(2),idecal(3),idecal(4),idecal(5)) 
279     &                               == Agrif_SpecialValue) then
280                           imin(iii) = imax(iii)
281                        endif 
282                     CASE (6)
283                        if (tempP%var%array6(idecal(1),
284     &            idecal(2),idecal(3),idecal(4),idecal(5),idecal(6)) 
285     &                               == Agrif_SpecialValue) then
286                           imin(iii) = imax(iii)
287                        endif                                                                                           
288                     END SELECT
289                  endif
290                  endif
291               endif           
292            enddo
293C
294            Res = 0.
295            Nbvals = 0
296C
297            SELECT CASE(nbdim)
298            CASE (1)
299!CDIR ALTCODE
300!CDIR SHORTLOOP
301               do ii = imin(1),imax(1)
302                    ValParent = parent%var%array1(ii)
303                    if ( ValParent .NE. Agrif_SpecialValue) then
304                        Res = Res + ValParent
305                        Nbvals = Nbvals + 1
306                    endif
307               enddo
308C
309            CASE (2)
310               do jj = imin(2),imax(2)
311!CDIR ALTCODE
312!CDIR SHORTLOOP
313               do ii = imin(1),imax(1)
314                    ValParent = parent%var%array2(ii,jj)
315                    if ( ValParent .NE. Agrif_SpecialValue) then
316                        Res = Res + ValParent
317                        Nbvals = Nbvals + 1
318                    endif
319               enddo 
320               enddo
321               
322            CASE (3)
323               do kk = imin(3),imax(3)
324               do jj = imin(2),imax(2)
325!CDIR ALTCODE
326!CDIR SHORTLOOP
327               do ii = imin(1),imax(1)
328                    ValParent = parent%var%array3(ii,jj,kk)
329                    if ( ValParent .NE. Agrif_SpecialValue) then
330                        Res = Res + ValParent
331                        Nbvals = Nbvals + 1
332                    endif
333                        enddo
334                  enddo 
335               enddo
336
337            CASE (4)
338               do ll = imin(4),imax(4)
339               do kk = imin(3),imax(3)
340               do jj = imin(2),imax(2)
341!CDIR ALTCODE
342!CDIR SHORTLOOP
343               do ii = imin(1),imax(1)
344                    ValParent = parent%var%array4(ii,jj,kk,ll)
345                    if ( ValParent .NE. Agrif_SpecialValue) then
346                        Res = Res + ValParent
347                        Nbvals = Nbvals + 1
348                    endif
349                              enddo
350                        enddo
351                  enddo 
352               enddo
353
354            CASE (5)
355               do mm = imin(5),imax(5)
356               do ll = imin(4),imax(4)
357               do kk = imin(3),imax(3)
358               do jj = imin(2),imax(2)
359!CDIR ALTCODE
360!CDIR SHORTLOOP
361               do ii = imin(1),imax(1)
362                    ValParent = parent%var%array5(ii,jj,kk,ll,mm)
363                    if ( ValParent .NE. Agrif_SpecialValue) then
364                        Res = Res + ValParent
365                        Nbvals = Nbvals + 1
366                    endif
367                                    enddo
368                              enddo
369                        enddo
370                  enddo 
371               enddo
372
373            CASE (6)
374               do nn = imin(6),imax(6)
375               do mm = imin(5),imax(5)
376               do ll = imin(4),imax(4)
377               do kk = imin(3),imax(3)
378               do jj = imin(2),imax(2)
379!CDIR ALTCODE
380!CDIR SHORTLOOP
381               do ii = imin(1),imax(1)
382                    ValParent = parent%var%array6(ii,jj,kk,ll,mm,nn)
383                    if ( ValParent .NE. Agrif_SpecialValue) then
384                        Res = Res + ValParent
385                        Nbvals = Nbvals + 1
386                    endif
387                                          enddo
388                                    enddo
389                              enddo
390                        enddo
391                  enddo 
392               enddo
393
394            END SELECT
395C
396C
397           
398            if (Nbvals.GT.0) then
399              if (firsttest) then
400                   firsttest = .FALSE.
401                   i=1
402                   cycle
403              endif
404            SELECT CASE(nbdim)
405            CASE (1)             
406              tempP%var%array1(indic(1)) 
407     &           = Res/Nbvals
408            CASE (2)
409              tempP%var%array2(indic(1),
410     &                            indic(2)) = Res/Nbvals
411            CASE (3)
412              tempP%var%array3(indic(1),
413     &                            indic(2),indic(3)) = Res/Nbvals
414            CASE (4)
415              tempP%var%array4(indic(1),
416     &                            indic(2),indic(3),indic(4))
417     &                = Res/Nbvals
418            CASE (5)
419              tempP%var%array5(indic(1),
420     &                            indic(2),indic(3),indic(4),
421     &                   indic(5)) = Res/Nbvals
422            CASE (6)
423              tempP%var%array6(indic(1),
424     &                            indic(2),indic(3),indic(4),
425     &                           indic(5),indic(6)) = Res/Nbvals
426            END SELECT
427              exit
428            else
429               if (firsttest) exit
430               i = i + 1                     
431            endif
432          enddo           
433C     
434      End Subroutine CalculNewValTempP
435C
436C
437      End Module Agrif_Mask     
438     
439      Subroutine CalculNewValTempP3D(indic,
440     &               tempP,parent,ppbtab,
441     &               ppetab,noraftab,MaxSearch,Agrif_SpecialValue)
442C
443CCC   Description:
444CCC   Subroutine called in the procedure Agrif_InterpnD to recalculate the value
445CCC   of the parent grid variable when this one is equal to Agrif_SpecialValue. 
446C
447C     Declarations:
448C
449       
450C
451C     Arrays arguments     
452      INTEGER, PARAMETER :: nbdim = 3
453      INTEGER,DIMENSION(nbdim) :: indic
454      LOGICAL,DIMENSION(nbdim) :: noraftab
455      INTEGER,DIMENSION(nbdim) :: ppbtab,ppetab
456      REAL :: Agrif_SpecialValue
457C
458C     Pointer argument
459      REAL,DIMENSION(ppbtab(1):ppetab(1),ppbtab(2):ppetab(2),
460     &    ppbtab(3):ppetab(3)) :: tempP, parent  ! Part of the parent grid used for
461                                      ! the interpolation of the child grid                     
462C
463C     Local scalar
464      INTEGER                  :: i,ii,iii,jj,kk,ll,mm,nn 
465      INTEGER,DIMENSION(nbdim) :: imin,imax,idecal
466      INTEGER                  :: Nbvals
467      REAL                     :: Res
468      REAL                     :: ValParent
469      INTEGER                  :: ValMax
470      INTEGER :: MaxSearch
471      LOGICAL :: Existunmasked
472C     
473C     Local arrays     
474C
475      ValMax = 1
476!CDIR NOVECTOR
477      do iii = 1 , nbdim
478         IF (.NOT.noraftab(iii)) THEN
479         ValMax = max(ValMax,ppetab(iii)-indic(iii))
480         ValMax = max(ValMax,indic(iii)-ppbtab(iii))
481         ENDIF
482      enddo
483C
484      Valmax = min(Valmax,MaxSearch)
485C
486!CDIR NOVECTOR
487      imin = indic
488!CDIR NOVECTOR
489      imax = indic
490     
491!CDIR NOVECTOR
492         idecal = indic 
493C
494         i = Valmax
495
496            do iii = 1 , nbdim
497               if (.NOT.noraftab(iii)) then
498                  imin(iii) = max(indic(iii) - i,ppbtab(iii))
499                  imax(iii) = min(indic(iii) + i,ppetab(iii))
500         
501                  if (indic(iii).GT.ppbtab(iii)) then
502
503                     idecal(iii) = idecal(iii)-1
504
505                        if (tempP(idecal(1),idecal(2),idecal(3)) 
506     &                               == Agrif_SpecialValue) then
507                           imin(iii) = imax(iii)
508                        endif
509
510                     idecal(iii) = idecal(iii)+1
511                  endif
512
513               endif           
514            enddo
515C
516            Existunmasked = .FALSE.
517C
518               do kk = imin(3),imax(3)
519               do jj = imin(2),imax(2)
520!CDIR NOVECTOR
521               do ii = imin(1),imax(1)
522                    if ( parent(ii,jj,kk) .NE. Agrif_SpecialValue) then
523                        Existunmasked = .TRUE.
524                       exit
525                    endif
526                        enddo
527                  enddo 
528               enddo
529C
530C
531          if (.Not.Existunmasked) return
532C
533         i = 1
534C
535         do While(i <= ValMax)
536C
537
538            do iii = 1 , nbdim
539               if (.NOT.noraftab(iii)) then
540                  imin(iii) = max(indic(iii) - i,ppbtab(iii))
541                  imax(iii) = min(indic(iii) + i,ppetab(iii))
542               endif           
543            enddo
544C
545            Res = 0.
546            Nbvals = 0
547C
548               do kk = imin(3),imax(3)
549               do jj = imin(2),imax(2)
550!CDIR NOVECTOR
551               do ii = imin(1),imax(1)
552                    ValParent = parent(ii,jj,kk)
553                    if ( ValParent .NE. Agrif_SpecialValue) then
554                        Res = Res + ValParent
555                        Nbvals = Nbvals + 1
556                    endif
557                        enddo
558                  enddo 
559               enddo
560C
561C
562           
563            if (Nbvals.GT.0) then
564              tempP(indic(1),indic(2),indic(3)) = Res/Nbvals
565              exit
566            else
567               i = i + 1                     
568            endif
569          enddo           
570C     
571      End Subroutine CalculNewValTempP3D       
572
573      Subroutine CalculNewValTempP4D(indic,
574     &               tempP,parent,ppbtab,
575     &               ppetab,noraftab,MaxSearch,Agrif_SpecialValue)
576C
577CCC   Description:
578CCC   Subroutine called in the procedure Agrif_InterpnD to recalculate the value
579CCC   of the parent grid variable when this one is equal to Agrif_SpecialValue. 
580C
581C     Declarations:
582C
583       
584C
585C     Arrays arguments     
586      INTEGER, PARAMETER :: nbdim = 4
587      INTEGER,DIMENSION(nbdim) :: indic
588      LOGICAL,DIMENSION(nbdim) :: noraftab
589      INTEGER,DIMENSION(nbdim) :: ppbtab,ppetab
590      INTEGER :: MaxSearch     
591      REAL :: Agrif_SpecialValue     
592C
593C     Pointer argument
594      REAL,DIMENSION(ppbtab(1):ppetab(1),ppbtab(2):ppetab(2),
595     &    ppbtab(3):ppetab(3),
596     &    ppbtab(4):ppetab(4)) :: tempP, parent  ! Part of the parent grid used for
597                                      ! the interpolation of the child grid                     
598C
599C     Local scalar
600      INTEGER                  :: i,ii,iii,jj,kk,ll,mm,nn 
601      INTEGER,DIMENSION(nbdim) :: imin,imax,idecal
602      INTEGER                  :: Nbvals
603      REAL                     :: Res
604      REAL                     :: ValParent
605      INTEGER                  :: ValMax
606      LOGICAL                  :: firsttest
607C     
608C     Local arrays     
609C
610      ValMax = 1
611      do iii = 1 , nbdim
612         IF (.NOT.noraftab(iii)) THEN
613         ValMax = max(ValMax,ppetab(iii)-indic(iii))
614         ValMax = max(ValMax,indic(iii)-ppbtab(iii))
615         ENDIF
616      enddo
617C
618      Valmax = min(Valmax,MaxSearch)
619C
620      imin = indic
621      imax = indic
622C
623         i = 1
624         firsttest = .TRUE.
625         idecal = indic 
626   
627C
628         do While(i <= ValMax)
629C
630         IF ((i == 1).AND.(firsttest)) i = Valmax
631
632            do iii = 1 , nbdim
633               if (.NOT.noraftab(iii)) then
634                  imin(iii) = max(indic(iii) - i,ppbtab(iii))
635                  imax(iii) = min(indic(iii) + i,ppetab(iii))
636                  if (firsttest) then               
637                  if (indic(iii).GT.ppbtab(iii)) then
638
639
640                     idecal(iii) = idecal(iii)-1
641           
642                     if (tempP(idecal(1),idecal(2),idecal(3),idecal(4))
643     &                               == Agrif_SpecialValue) then
644                           imin(iii) = imax(iii)
645                     endif
646           
647           idecal(iii) = idecal(iii)+1
648                  endif
649                  endif
650               endif           
651            enddo
652C
653            Res = 0.
654            Nbvals = 0
655C
656               do ll = imin(4),imax(4)
657               do kk = imin(3),imax(3)
658               do jj = imin(2),imax(2)
659               do ii = imin(1),imax(1)
660                    ValParent = parent(ii,jj,kk,ll)
661                    if ( ValParent .NE. Agrif_SpecialValue) then
662                        Res = Res + ValParent
663                        Nbvals = Nbvals + 1
664                    endif
665                              enddo
666                        enddo
667                  enddo 
668               enddo
669C
670C
671           
672            if (Nbvals.GT.0) then
673              if (firsttest) then
674                   firsttest = .FALSE.
675                   i=1
676                   cycle
677              endif
678         
679              tempP(indic(1),indic(2),indic(3),indic(4)) = Res/Nbvals
680       
681              exit
682            else
683               if (firsttest) exit
684               i = i + 1                     
685            endif
686          enddo           
687C     
688      End Subroutine CalculNewValTempP4D       
Note: See TracBrowser for help on using the repository browser.