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 trunk/AGRIF/AGRIF_FILES – NEMO

source: trunk/AGRIF/AGRIF_FILES/modmask.F @ 396

Last change on this file since 396 was 396, checked in by opalod, 18 years ago

Initial revision

  • Property svn:eol-style set to native
  • Property svn:keywords set to Author Date Id Revision
File size: 12.0 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
37      Integer, Parameter :: MaxSearch = 5
38C
39      CONTAINS
40C     Define procedures contained in this module
41C
42C     **************************************************************************
43C     Subroutine Agrif_CheckMasknD
44C     **************************************************************************
45C       
46      Subroutine Agrif_CheckMasknD(tempP,parent,pbtab,petab,ppbtab,
47     &               ppetab,noraftab,nbdim)
48C
49CCC   Description:
50CCC   Subroutine called in the procedure Agrif_InterpnD to recalculate the value
51CCC   of the parent grid variable when this one is equal to Agrif_SpecialValue. 
52C
53C     Declarations:
54C
55       
56C
57C     Arrays arguments     
58      INTEGER :: nbdim
59      INTEGER,DIMENSION(nbdim) :: pbtab  ! Limits of the parent grid used 
60      INTEGER,DIMENSION(nbdim) :: petab  ! interpolation of the child grid
61      LOGICAL,DIMENSION(nbdim) :: noraftab
62      INTEGER,DIMENSION(nbdim) :: ppbtab,ppetab
63C
64C     Pointer argument
65      TYPE(AGRIF_PVARIABLE) :: tempP  ! Part of the parent grid used for
66                                      ! the interpolation of the child grid                     
67C     Data TYPE argument                                   
68      TYPE(AGRIF_PVARIABLE) :: parent      ! The parent grid
69C
70C     Local scalar
71      INTEGER                   :: i0,j0,k0,l0,m0,n0
72C     
73C     Local arrays
74C
75C     
76      SELECT CASE (nbdim)
77      CASE (1)
78         do i0 = pbtab(1),petab(1)
79         IF (tempP%var%array1(i0)
80     &                        == Agrif_SpecialValue) Then
81            Call CalculNewValTempP((/i0/),
82     &                             tempP,parent,
83     &                             ppbtab,ppetab,
84     &                             noraftab,nbdim)
85         ENDIF
86         enddo
87      CASE (2)
88         do j0 = pbtab(2),petab(2)
89         do i0 = pbtab(1),petab(1)
90         IF (tempP%var%array2(i0,j0)
91     &                        == Agrif_SpecialValue) Then
92            Call CalculNewValTempP((/i0,j0/),
93     &                             tempP,parent,
94     &                             ppbtab,ppetab,
95     &                             noraftab,nbdim)
96         ENDIF
97         enddo 
98         enddo
99      CASE (3)
100         do k0 = pbtab(3),petab(3)
101         do j0 = pbtab(2),petab(2)
102         do i0 = pbtab(1),petab(1)
103         IF (tempP%var%array3(i0,j0,k0)
104     &                        == Agrif_SpecialValue) Then
105            Call CalculNewValTempP((/i0,j0,k0/),
106     &                             tempP,parent,
107     &                             ppbtab,ppetab,
108     &                             noraftab,nbdim)
109         ENDIF
110         enddo
111         enddo 
112         enddo
113      CASE (4)
114         do l0 = pbtab(4),petab(4)
115         do k0 = pbtab(3),petab(3)
116         do j0 = pbtab(2),petab(2)
117         do i0 = pbtab(1),petab(1)
118         IF (tempP%var%array4(i0,j0,k0,l0) 
119     &                        == Agrif_SpecialValue) Then
120            Call CalculNewValTempP((/i0,j0,k0,l0/),
121     &                             tempP,parent,
122     &                             ppbtab,ppetab,
123     &                             noraftab,nbdim)
124         ENDIF
125         enddo
126         enddo
127         enddo 
128         enddo
129      CASE (5)
130         do m0 = pbtab(5),petab(5)
131         do l0 = pbtab(4),petab(4)
132         do k0 = pbtab(3),petab(3)
133         do j0 = pbtab(2),petab(2)
134         do i0 = pbtab(1),petab(1)
135         IF (tempP%var%array5(i0,j0,k0,l0,m0) 
136     &                       == Agrif_SpecialValue) Then
137            Call CalculNewValTempP((/i0,j0,k0,l0,m0/),
138     &                             tempP,parent,
139     &                             ppbtab,ppetab,
140     &                             noraftab,nbdim)
141         ENDIF
142         enddo
143         enddo
144         enddo 
145         enddo
146         enddo
147      CASE (6)
148         do n0 = pbtab(6),petab(6)
149         do m0 = pbtab(5),petab(5)
150         do l0 = pbtab(4),petab(4)
151         do k0 = pbtab(3),petab(3)
152         do j0 = pbtab(2),petab(2)
153         do i0 = pbtab(1),petab(1)
154         IF (tempP%var%array6(i0,j0,k0,l0,m0,n0) 
155     &                       == Agrif_SpecialValue) Then
156            Call CalculNewValTempP((/i0,j0,k0,l0,m0,n0/),
157     &                             tempP,parent,
158     &                             ppbtab,ppetab,
159     &                             noraftab,nbdim)
160         ENDIF
161         enddo
162         enddo
163         enddo
164         enddo 
165         enddo
166         enddo
167      END SELECT
168C
169C     
170C     
171      End Subroutine Agrif_CheckMasknD
172C
173C
174C     **************************************************************************
175C     Subroutine CalculNewValTempP
176C     **************************************************************************
177C       
178      Subroutine CalculNewValTempP(indic,
179     &               tempP,parent,ppbtab,
180     &               ppetab,noraftab,nbdim)
181C
182CCC   Description:
183CCC   Subroutine called in the procedure Agrif_InterpnD to recalculate the value
184CCC   of the parent grid variable when this one is equal to Agrif_SpecialValue. 
185C
186C     Declarations:
187C
188       
189C
190C     Arrays arguments     
191      INTEGER :: nbdim
192      INTEGER,DIMENSION(nbdim) :: indic
193      LOGICAL,DIMENSION(nbdim) :: noraftab
194      INTEGER,DIMENSION(nbdim) :: ppbtab,ppetab
195C
196C     Pointer argument
197      TYPE(AGRIF_PVARIABLE) :: tempP  ! Part of the parent grid used for
198                                      ! the interpolation of the child grid                     
199C     Data TYPE argument                                   
200      TYPE(AGRIF_PVARIABLE) :: parent      ! The parent grid
201C
202C     Local scalar
203      INTEGER                  :: i,ii,iii,jj,kk,ll,mm,nn 
204      INTEGER,DIMENSION(nbdim) :: imin,imax
205      INTEGER                  :: Nbvals
206      REAL                     :: Res
207      REAL                     :: ValParent
208      INTEGER                  :: ValMax
209      LOGICAL                  :: firsttest
210C     
211C     Local arrays     
212C
213      ValMax = 1
214      do iii = 1 , nbdim
215         IF (.NOT.noraftab(iii)) THEN
216         ValMax = max(ValMax,ppetab(iii)-indic(iii))
217         ValMax = max(ValMax,indic(iii)-ppbtab(iii))
218         ENDIF
219      enddo
220C
221      Valmax = min(Valmax,MaxSearch)
222C
223      imin = indic
224      imax = indic
225C
226         i = 1
227         firsttest = .TRUE.
228C
229         do While(i <= ValMax)
230C
231         IF ((i == 1).AND.(firsttest)) i = Valmax
232
233            do iii = 1 , nbdim
234               if (.NOT.noraftab(iii)) then
235                  imin(iii) = max(indic(iii) - i,ppbtab(iii))
236                  imax(iii) = min(indic(iii) + i,ppetab(iii))
237               endif           
238            enddo
239C
240            Res = 0.
241            Nbvals = 0
242C
243            if ( nbdim .EQ. 1 ) then
244               do ii = imin(1),imax(1)
245                    ValParent = parent%var%array1(ii)
246                    if ( ValParent .NE. Agrif_SpecialValue) then
247                        Res = Res + ValParent
248                        Nbvals = Nbvals + 1
249                    endif
250               enddo
251            endif
252C
253            if ( nbdim .EQ. 2 ) then
254               do jj = imin(2),imax(2)
255               do ii = imin(1),imax(1)
256                    ValParent = parent%var%array2(ii,jj)
257                    if ( ValParent .NE. Agrif_SpecialValue) then
258                        Res = Res + ValParent
259                        Nbvals = Nbvals + 1
260                    endif
261               enddo 
262               enddo
263            endif
264C
265            if ( nbdim .EQ. 3 ) then
266               do kk = imin(3),imax(3)
267               do jj = imin(2),imax(2)
268               do ii = imin(1),imax(1)
269                    ValParent = parent%var%array3(ii,jj,kk)
270                    if ( ValParent .NE. Agrif_SpecialValue) then
271                        Res = Res + ValParent
272                        Nbvals = Nbvals + 1
273                    endif
274                        enddo
275                  enddo 
276               enddo
277            endif
278C
279            if ( nbdim .EQ. 4 ) then
280               do ll = imin(4),imax(4)
281               do kk = imin(3),imax(3)
282               do jj = imin(2),imax(2)
283               do ii = imin(1),imax(1)
284                    ValParent = parent%var%array4(ii,jj,kk,ll)
285                    if ( ValParent .NE. Agrif_SpecialValue) then
286                        Res = Res + ValParent
287                        Nbvals = Nbvals + 1
288                    endif
289                              enddo
290                        enddo
291                  enddo 
292               enddo
293            endif
294C
295            if ( nbdim .EQ. 5 ) then
296               do mm = imin(5),imax(5)
297               do ll = imin(4),imax(4)
298               do kk = imin(3),imax(3)
299               do jj = imin(2),imax(2)
300               do ii = imin(1),imax(1)
301                    ValParent = parent%var%array5(ii,jj,kk,ll,mm)
302                    if ( ValParent .NE. Agrif_SpecialValue) then
303                        Res = Res + ValParent
304                        Nbvals = Nbvals + 1
305                    endif
306                                    enddo
307                              enddo
308                        enddo
309                  enddo 
310               enddo
311            endif
312C
313            if ( nbdim .EQ. 6 ) then
314               do nn = imin(6),imax(6)
315               do mm = imin(5),imax(5)
316               do ll = imin(4),imax(4)
317               do kk = imin(3),imax(3)
318               do jj = imin(2),imax(2)
319               do ii = imin(1),imax(1)
320                    ValParent = parent%var%array6(ii,jj,kk,ll,mm,nn)
321                    if ( ValParent .NE. Agrif_SpecialValue) then
322                        Res = Res + ValParent
323                        Nbvals = Nbvals + 1
324                    endif
325                                          enddo
326                                    enddo
327                              enddo
328                        enddo
329                  enddo 
330               enddo
331            endif
332C
333C
334           
335            if (Nbvals.GT.0) then
336              if (firsttest) then
337                   firsttest = .FALSE.
338                   cycle
339              endif
340              if ( nbdim .EQ. 1 ) tempP%var%array1(indic(1)) 
341     &           = Res/Nbvals
342              if ( nbdim .EQ. 2 ) tempP%var%array2(indic(1),
343     &                            indic(2)) = Res/Nbvals
344              if ( nbdim .EQ. 3 ) tempP%var%array3(indic(1),
345     &                            indic(2),indic(3)) = Res/Nbvals
346              if ( nbdim .EQ. 4 ) tempP%var%array4(indic(1),
347     &                            indic(2),indic(3),indic(4))
348     &                = Res/Nbvals
349              if ( nbdim .EQ. 5 ) tempP%var%array5(indic(1),
350     &                            indic(2),indic(3),indic(4),
351     &                   indic(5)) = Res/Nbvals
352              if ( nbdim .EQ. 6 ) tempP%var%array6(indic(1),
353     &                            indic(2),indic(3),indic(4),
354     &                           indic(5),indic(6)) = Res/Nbvals
355              exit
356            else
357               if (firsttest) exit
358               i = i + 1                     
359            endif
360          enddo           
361C     
362      End Subroutine CalculNewValTempP
363C
364C
365      End Module Agrif_Mask       
Note: See TracBrowser for help on using the repository browser.