source: codes/icosagcm/devel/src/base/field.f90 @ 595

Last change on this file since 595 was 533, checked in by dubos, 7 years ago

devel : reorganization of source tree

File size: 14.7 KB
Line 
1MODULE field_mod
2  USE genmod
3  IMPLICIT NONE
4 
5  INTEGER,PARAMETER :: field_T=1
6  INTEGER,PARAMETER :: field_U=2
7  INTEGER,PARAMETER :: field_Z=3
8
9  INTEGER,PARAMETER :: type_integer=1
10  INTEGER,PARAMETER :: type_real=2
11  INTEGER,PARAMETER :: type_logical=3
12   
13  TYPE t_field
14    CHARACTER(30)      :: name
15    REAL(rstd),POINTER :: rval2d(:)
16    REAL(rstd),POINTER :: rval3d(:,:)
17    REAL(rstd),POINTER :: rval4d(:,:,:)
18
19    INTEGER,POINTER :: ival2d(:)
20    INTEGER,POINTER :: ival3d(:,:)
21    INTEGER,POINTER :: ival4d(:,:,:)
22   
23    LOGICAL,POINTER :: lval2d(:)
24    LOGICAL,POINTER :: lval3d(:,:)
25    LOGICAL,POINTER :: lval4d(:,:,:)
26
27    INTEGER :: ndim
28    INTEGER :: field_type
29    INTEGER :: data_type 
30    INTEGER :: dim3
31    INTEGER :: dim4
32  END TYPE t_field   
33
34  INTERFACE get_val
35    MODULE PROCEDURE getval_r2d,getval_r3d,getval_r4d, &
36                     getval_i2d,getval_i3d,getval_i4d, &
37                     getval_l2d,getval_l3d,getval_l4d
38  END INTERFACE
39                   
40  INTERFACE ASSIGNMENT(=)
41    MODULE PROCEDURE getval_r2d,getval_r3d,getval_r4d, &
42                     getval_i2d,getval_i3d,getval_i4d, &
43                     getval_l2d,getval_l3d,getval_l4d 
44  END INTERFACE
45
46  PRIVATE :: allocate_field_
47
48CONTAINS
49
50  SUBROUTINE allocate_field(field,field_type,data_type,dim1,dim2,name)
51  USE domain_mod
52  USE omp_para
53    TYPE(t_field),POINTER :: field(:)
54    INTEGER,INTENT(IN) :: field_type
55    INTEGER,INTENT(IN) :: data_type
56    INTEGER,OPTIONAL :: dim1,dim2
57    CHARACTER(*), OPTIONAL :: name
58!$OMP BARRIER
59!$OMP MASTER
60    ALLOCATE(field(ndomain))
61!$OMP END MASTER
62!$OMP BARRIER
63    CALL allocate_field_(field,field_type,data_type,dim1,dim2,name)
64  END SUBROUTINE allocate_field
65
66  SUBROUTINE allocate_fields(nfield,field,field_type,data_type,dim1,dim2,name)
67  USE domain_mod
68  USE omp_para
69    INTEGER,INTENT(IN) :: nfield
70    TYPE(t_field),POINTER :: field(:,:)
71    INTEGER,INTENT(IN) :: field_type
72    INTEGER,INTENT(IN) :: data_type
73    INTEGER,OPTIONAL :: dim1,dim2
74    CHARACTER(*), OPTIONAL :: name
75    INTEGER :: i
76!$OMP BARRIER
77!$OMP MASTER
78    ALLOCATE(field(ndomain,nfield))
79!$OMP END MASTER
80!$OMP BARRIER
81    DO i=1,nfield
82       CALL allocate_field_(field(:,i),field_type,data_type,dim1,dim2,name)
83    END DO
84  END SUBROUTINE allocate_fields
85
86  SUBROUTINE allocate_field_(field,field_type,data_type,dim1,dim2,name)
87  USE domain_mod
88  USE omp_para
89  IMPLICIT NONE
90    TYPE(t_field) :: field(:)
91    INTEGER,INTENT(IN) :: field_type
92    INTEGER,INTENT(IN) :: data_type
93    INTEGER,OPTIONAL :: dim1,dim2
94    CHARACTER(*), OPTIONAL :: name
95    INTEGER :: ind
96    INTEGER :: ii_size,jj_size
97
98    DO ind=1,ndomain
99      IF (.NOT. assigned_domain(ind)  .OR. .NOT. is_omp_level_master) CYCLE
100
101      IF(PRESENT(name)) THEN
102         field(ind)%name = name
103      ELSE
104         field(ind)%name = '(undefined)'
105      END IF
106
107      IF (PRESENT(dim2)) THEN
108        field(ind)%ndim=4 
109        field(ind)%dim4=dim2 
110        field(ind)%dim3=dim1
111      ELSE IF (PRESENT(dim1)) THEN
112        field(ind)%ndim=3
113        field(ind)%dim3=dim1
114      ELSE
115        field(ind)%ndim=2
116      ENDIF
117   
118   
119      field(ind)%data_type=data_type
120      field(ind)%field_type=field_type
121   
122      IF (field_type==field_T) THEN
123        jj_size=domain(ind)%jjm
124      ELSE IF (field_type==field_U) THEN
125        jj_size=3*domain(ind)%jjm
126      ELSE IF (field_type==field_Z) THEN
127        jj_size=2*domain(ind)%jjm
128      ENDIF
129     
130      ii_size=domain(ind)%iim
131       
132      IF (field(ind)%ndim==4) THEN
133        IF (data_type==type_integer) ALLOCATE(field(ind)%ival4d(ii_size*jj_size,dim1,dim2))
134        IF (data_type==type_real)    ALLOCATE(field(ind)%rval4d(ii_size*jj_size,dim1,dim2))
135        IF (data_type==type_logical) ALLOCATE(field(ind)%lval4d(ii_size*jj_size,dim1,dim2))
136      ELSE IF (field(ind)%ndim==3) THEN
137        IF (data_type==type_integer) ALLOCATE(field(ind)%ival3d(ii_size*jj_size,dim1))
138        IF (data_type==type_real)    ALLOCATE(field(ind)%rval3d(ii_size*jj_size,dim1))
139        IF (data_type==type_logical) ALLOCATE(field(ind)%lval3d(ii_size*jj_size,dim1))
140      ELSE IF (field(ind)%ndim==2) THEN
141        IF (data_type==type_integer) ALLOCATE(field(ind)%ival2d(ii_size*jj_size))
142        IF (data_type==type_real)    ALLOCATE(field(ind)%rval2d(ii_size*jj_size))
143        IF (data_type==type_logical) ALLOCATE(field(ind)%lval2d(ii_size*jj_size))
144      ENDIF
145     
146   ENDDO
147!$OMP BARRIER
148   
149 END SUBROUTINE allocate_field_
150
151  SUBROUTINE allocate_field_glo(field,field_type,data_type,dim1,dim2,name)
152  USE domain_mod
153  IMPLICIT NONE
154    TYPE(t_field),POINTER :: field(:)
155    INTEGER,INTENT(IN) :: field_type
156    INTEGER,INTENT(IN) :: data_type
157    INTEGER,OPTIONAL :: dim1,dim2
158    CHARACTER(*), OPTIONAL :: name
159    INTEGER :: ind
160    INTEGER :: ii_size,jj_size
161
162    ALLOCATE(field(ndomain_glo))   
163
164    DO ind=1,ndomain_glo
165 
166      IF (PRESENT(dim2)) THEN
167        field(ind)%ndim=4 
168        field(ind)%dim4=dim2 
169        field(ind)%dim3=dim1 
170      ELSE IF (PRESENT(dim1)) THEN
171        field(ind)%ndim=3
172        field(ind)%dim3=dim1 
173      ELSE
174        field(ind)%ndim=2
175      ENDIF
176   
177      IF(PRESENT(name)) THEN
178         field(ind)%name = name
179      ELSE
180         field(ind)%name = '(undefined)'
181      END IF
182   
183      field(ind)%data_type=data_type
184      field(ind)%field_type=field_type
185   
186      IF (field_type==field_T) THEN
187        jj_size=domain_glo(ind)%jjm
188      ELSE IF (field_type==field_U) THEN
189        jj_size=3*domain_glo(ind)%jjm
190      ELSE IF (field_type==field_Z) THEN
191        jj_size=2*domain_glo(ind)%jjm
192      ENDIF
193     
194      ii_size=domain_glo(ind)%iim
195       
196      IF (field(ind)%ndim==4) THEN
197        IF (data_type==type_integer) ALLOCATE(field(ind)%ival4d(ii_size*jj_size,dim1,dim2))
198        IF (data_type==type_real)    ALLOCATE(field(ind)%rval4d(ii_size*jj_size,dim1,dim2))
199        IF (data_type==type_logical) ALLOCATE(field(ind)%lval4d(ii_size*jj_size,dim1,dim2))
200      ELSE IF (field(ind)%ndim==3) THEN
201        IF (data_type==type_integer) ALLOCATE(field(ind)%ival3d(ii_size*jj_size,dim1))
202        IF (data_type==type_real)    ALLOCATE(field(ind)%rval3d(ii_size*jj_size,dim1))
203        IF (data_type==type_logical) ALLOCATE(field(ind)%lval3d(ii_size*jj_size,dim1))
204      ELSE IF (field(ind)%ndim==2) THEN
205        IF (data_type==type_integer) ALLOCATE(field(ind)%ival2d(ii_size*jj_size))
206        IF (data_type==type_real)    ALLOCATE(field(ind)%rval2d(ii_size*jj_size))
207        IF (data_type==type_logical) ALLOCATE(field(ind)%lval2d(ii_size*jj_size))
208      ENDIF
209     
210   ENDDO
211 
212  END SUBROUTINE allocate_field_glo
213
214  SUBROUTINE deallocate_field(field)
215    USE domain_mod
216    USE omp_para
217    IMPLICIT NONE
218    TYPE(t_field),POINTER :: field(:)
219    !$OMP BARRIER
220    CALL deallocate_field_(field)
221    !$OMP BARRIER
222    !$OMP MASTER
223    DEALLOCATE(field)
224    !$OMP END MASTER
225    !$OMP BARRIER
226  END SUBROUTINE deallocate_field
227 
228  SUBROUTINE deallocate_fields(field)
229    USE domain_mod
230    USE omp_para
231    IMPLICIT NONE
232    TYPE(t_field),POINTER :: field(:,:)
233    INTEGER :: i
234    !$OMP BARRIER
235    DO i=1,SIZE(field,2)
236       CALL deallocate_field_(field(:,i))
237    END DO
238    !$OMP BARRIER
239    !$OMP MASTER
240    DEALLOCATE(field)
241    !$OMP END MASTER
242    !$OMP BARRIER
243  END SUBROUTINE deallocate_fields
244
245  SUBROUTINE deallocate_field_(field)
246  USE domain_mod
247  USE omp_para
248  IMPLICIT NONE
249    TYPE(t_field) :: field(:)
250    INTEGER :: data_type
251    INTEGER :: ind
252    DO ind=1,ndomain
253      IF (.NOT. assigned_domain(ind)  .OR. .NOT. is_omp_level_master) CYCLE
254
255      data_type=field(ind)%data_type
256       
257      IF (field(ind)%ndim==4) THEN
258        IF (data_type==type_integer) DEALLOCATE(field(ind)%ival4d)
259        IF (data_type==type_real)    DEALLOCATE(field(ind)%rval4d)
260        IF (data_type==type_logical) DEALLOCATE(field(ind)%lval4d)
261      ELSE IF (field(ind)%ndim==3) THEN
262        IF (data_type==type_integer) DEALLOCATE(field(ind)%ival3d)
263        IF (data_type==type_real)    DEALLOCATE(field(ind)%rval3d)
264        IF (data_type==type_logical) DEALLOCATE(field(ind)%lval3d)
265      ELSE IF (field(ind)%ndim==2) THEN
266        IF (data_type==type_integer) DEALLOCATE(field(ind)%ival2d)
267        IF (data_type==type_real)    DEALLOCATE(field(ind)%rval2d)
268        IF (data_type==type_logical) DEALLOCATE(field(ind)%lval2d)
269      ENDIF
270     
271   ENDDO
272  END SUBROUTINE deallocate_field_
273
274  SUBROUTINE deallocate_field_glo(field)
275  USE domain_mod
276  IMPLICIT NONE
277    TYPE(t_field),POINTER :: field(:)
278    INTEGER :: data_type
279    INTEGER :: ind
280
281    DO ind=1,ndomain_glo
282
283      data_type=field(ind)%data_type
284       
285      IF (field(ind)%ndim==4) THEN
286        IF (data_type==type_integer) DEALLOCATE(field(ind)%ival4d)
287        IF (data_type==type_real)    DEALLOCATE(field(ind)%rval4d)
288        IF (data_type==type_logical) DEALLOCATE(field(ind)%lval4d)
289      ELSE IF (field(ind)%ndim==3) THEN
290        IF (data_type==type_integer) DEALLOCATE(field(ind)%ival3d)
291        IF (data_type==type_real)    DEALLOCATE(field(ind)%rval3d)
292        IF (data_type==type_logical) DEALLOCATE(field(ind)%lval3d)
293      ELSE IF (field(ind)%ndim==2) THEN
294        IF (data_type==type_integer) DEALLOCATE(field(ind)%ival2d)
295        IF (data_type==type_real)    DEALLOCATE(field(ind)%rval2d)
296        IF (data_type==type_logical) DEALLOCATE(field(ind)%lval2d)
297      ENDIF
298     
299   ENDDO
300   DEALLOCATE(field)
301       
302  END SUBROUTINE deallocate_field_glo
303   
304  SUBROUTINE extract_slice(field_in, field_out, l) 
305  USE domain_mod
306  USE omp_para
307  IMPLICIT NONE 
308    TYPE(t_field) :: field_in(:)
309    TYPE(t_field) :: field_out(:)
310    INTEGER,INTENT(IN) :: l
311   
312    INTEGER :: ind
313    INTEGER :: data_type
314
315!$OMP BARRIER
316    DO ind=1,ndomain
317      data_type=field_in(ind)%data_type
318      IF (.NOT. assigned_domain(ind)  .OR. .NOT. is_omp_level_master) CYCLE
319     
320      IF (field_in(ind)%ndim==3 .AND. field_out(ind)%ndim==2) THEN 
321        IF (data_type==type_integer)  field_out(ind)%ival2d=field_in(ind)%ival3d(:,l)
322        IF (data_type==type_real)     field_out(ind)%rval2d=field_in(ind)%rval3d(:,l)
323        IF (data_type==type_logical)  field_out(ind)%lval2d=field_in(ind)%lval3d(:,l)
324      ELSE IF  (field_in(ind)%ndim==4 .AND. field_out(ind)%ndim==3) THEN
325        IF (data_type==type_integer)  field_out(ind)%ival3d=field_in(ind)%ival4d(:,:,l)
326        IF (data_type==type_real)     field_out(ind)%rval3d=field_in(ind)%rval4d(:,:,l)
327        IF (data_type==type_logical)  field_out(ind)%lval3d=field_in(ind)%lval4d(:,:,l)
328      ELSE
329        PRINT *, 'extract_slice : cannot extract slice, dimension incompatible'
330        STOP       
331      ENDIF
332   ENDDO 
333!$OMP BARRIER   
334  END  SUBROUTINE extract_slice 
335 
336 
337  SUBROUTINE insert_slice(field_in, field_out, l) 
338  USE domain_mod
339  USE omp_para
340  IMPLICIT NONE 
341    TYPE(t_field) :: field_in(:)
342    TYPE(t_field) :: field_out(:)
343    INTEGER,INTENT(IN) :: l
344   
345    INTEGER :: ind
346    INTEGER :: data_type
347
348!$OMP BARRIER
349    DO ind=1,ndomain
350      data_type=field_in(ind)%data_type
351      IF (.NOT. assigned_domain(ind)  .OR. .NOT. is_omp_level_master) CYCLE
352     
353      IF (field_in(ind)%ndim==2 .AND. field_out(ind)%ndim==3) THEN 
354        IF (data_type==type_integer)  field_out(ind)%ival3d(:,l)=field_in(ind)%ival2d(:)
355        IF (data_type==type_real)     field_out(ind)%rval3d(:,l)=field_in(ind)%rval2d(:)
356        IF (data_type==type_logical)  field_out(ind)%lval3d(:,l)=field_in(ind)%lval2d(:)
357      ELSE IF  (field_in(ind)%ndim==4 .AND. field_out(ind)%ndim==3) THEN
358        IF (data_type==type_integer)  field_out(ind)%ival4d(:,:,l)=field_out(ind)%ival3d(:,:)
359        IF (data_type==type_real)     field_out(ind)%rval4d(:,:,l)=field_out(ind)%rval3d(:,:)
360        IF (data_type==type_logical)  field_out(ind)%lval4d(:,:,l)=field_out(ind)%lval3d(:,:)
361      ELSE
362        PRINT *, 'extract_slice : cannot insert slice, dimension incompatible'
363        STOP       
364      ENDIF
365   ENDDO 
366!$OMP BARRIER   
367 
368  END SUBROUTINE insert_slice
369   
370  SUBROUTINE getval_r2d(field_pt,field)
371  IMPLICIT NONE 
372    REAL(rstd),POINTER,INTENT(INOUT) :: field_pt(:)
373    TYPE(t_field),INTENT(IN) :: field
374   
375    IF (field%ndim/=2 .OR. field%data_type/=type_real) THEN
376       PRINT *, 'get_val_r2d : bad pointer assignment with ' // TRIM(field%name) 
377       STOP
378    END IF
379    field_pt=>field%rval2d
380  END SUBROUTINE  getval_r2d
381
382  SUBROUTINE getval_r3d(field_pt,field)
383  IMPLICIT NONE 
384    REAL(rstd),POINTER,INTENT(INOUT) :: field_pt(:,:)
385    TYPE(t_field),INTENT(IN) :: field
386   
387    IF (field%ndim/=3 .OR. field%data_type/=type_real) THEN
388       PRINT *, 'get_val_r3d : bad pointer assignment with ' // TRIM(field%name) 
389       STOP
390!       CALL ABORT
391    END IF
392    field_pt=>field%rval3d
393  END SUBROUTINE  getval_r3d
394
395  SUBROUTINE getval_r4d(field_pt,field)
396  IMPLICIT NONE 
397    REAL(rstd),POINTER,INTENT(INOUT) :: field_pt(:,:,:)
398    TYPE(t_field),INTENT(IN) :: field
399   
400    IF (field%ndim/=4 .OR. field%data_type/=type_real) THEN
401       PRINT *, 'get_val_r4d : bad pointer assignment with ' // TRIM(field%name)
402       STOP
403    END IF
404    field_pt=>field%rval4d
405  END SUBROUTINE  getval_r4d 
406
407 
408  SUBROUTINE getval_i2d(field_pt,field)
409  IMPLICIT NONE 
410    INTEGER,POINTER,INTENT(INOUT) :: field_pt(:)
411    TYPE(t_field),INTENT(IN) :: field
412   
413    IF (field%ndim/=2 .OR. field%data_type/=type_integer) STOP 'get_val_i2d : bad pointer assignment'       
414    field_pt=>field%ival2d
415  END SUBROUTINE  getval_i2d
416
417  SUBROUTINE getval_i3d(field_pt,field)
418  IMPLICIT NONE 
419    INTEGER,POINTER,INTENT(INOUT) :: field_pt(:,:)
420    TYPE(t_field),INTENT(IN) :: field
421   
422    IF (field%ndim/=3 .OR. field%data_type/=type_integer) STOP 'get_val_i3d : bad pointer assignment'       
423    field_pt=>field%ival3d
424  END SUBROUTINE  getval_i3d
425
426  SUBROUTINE getval_i4d(field_pt,field)
427  IMPLICIT NONE 
428    INTEGER,POINTER,INTENT(INOUT) :: field_pt(:,:,:)
429    TYPE(t_field),INTENT(IN) :: field
430   
431    IF (field%ndim/=4 .OR. field%data_type/=type_integer) STOP 'get_val_i4d : bad pointer assignment'       
432    field_pt=>field%ival4d
433  END SUBROUTINE  getval_i4d
434
435  SUBROUTINE getval_l2d(field_pt,field)
436  IMPLICIT NONE 
437    LOGICAL,POINTER,INTENT(INOUT) :: field_pt(:)
438    TYPE(t_field),INTENT(IN) :: field
439   
440    IF (field%ndim/=2 .OR. field%data_type/=type_logical) STOP 'get_val_l2d : bad pointer assignment'       
441    field_pt=>field%lval2d
442  END SUBROUTINE  getval_l2d
443
444  SUBROUTINE getval_l3d(field_pt,field)
445  IMPLICIT NONE 
446    LOGICAL,POINTER,INTENT(INOUT) :: field_pt(:,:)
447    TYPE(t_field),INTENT(IN) :: field
448   
449    IF (field%ndim/=3 .OR. field%data_type/=type_logical) STOP 'get_val_l3d : bad pointer assignment'       
450    field_pt=>field%lval3d
451  END SUBROUTINE  getval_l3d
452
453  SUBROUTINE getval_l4d(field_pt,field)
454  IMPLICIT NONE 
455    LOGICAL,POINTER,INTENT(INOUT) :: field_pt(:,:,:)
456    TYPE(t_field),INTENT(IN) :: field
457   
458    IF (field%ndim/=4 .OR. field%data_type/=type_logical) STOP 'get_val_l4d : bad pointer assignment'       
459    field_pt=>field%lval4d
460  END SUBROUTINE  getval_l4d   
461
462END MODULE field_mod   
Note: See TracBrowser for help on using the repository browser.