/[lmdze]/trunk/bibio/writefield.f
ViewVC logotype

Annotation of /trunk/bibio/writefield.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 113 - (hide annotations)
Thu Sep 18 19:56:46 2014 UTC (9 years, 8 months ago) by guez
File size: 2245 byte(s)
Moved the call to read_serre out of conf_gcm so that it can be called
only in the program ce0l, not in gcm. In gcm, variables of module
serre are read from start file. Added reading of dzoomx, dzoomy, taux,
tauy from start file, in dynetat0. Those variables were written by
dynredem0 but not read.

Removed possibility fxyhypb = false, because the geometric part of the
program is such a mess. Could then remove variables transx, transy,
alphax, alphay, pxo, pyo of module serre.

Bug fix in tau2alpha: missing save attributes. The first call to
tau2alpha needs to compute dxdyu and dxdyv regardless of value of
argument type, because they will be needed for subsequent calls to
tau2alpha with various values of argument type.

1 guez 108 module WriteField_m
2    
3 guez 110 use CreateNewField_m, only: CreateNewField, ncid, nbfield, MaxWriteField
4 guez 109 use GetFieldIndex_m, only: GetFieldIndex
5     USE netcdf95, ONLY: nf95_put_var
6 guez 108
7     implicit none
8    
9 guez 110 integer, save:: Record(MaxWriteField)
10    
11 guez 108 interface WriteField
12     module procedure WriteField3d,WriteField2d,WriteField1d
13     end interface WriteField
14    
15     private
16     public WriteField
17    
18     contains
19    
20     subroutine WriteField1d(name,Field)
21    
22 guez 109 character(len=*), intent(in):: name
23     real, intent(in):: Field(:)
24 guez 108
25 guez 109 ! Local:
26     integer index
27    
28     !-------------------------------------------
29    
30     Index=GetFieldIndex(name)
31    
32     if (Index==-1) then
33     call CreateNewField(name, shape(field))
34     Index=nbfield
35 guez 110 Record(Index) = 1
36 guez 109 else
37     Record(Index)=Record(Index)+1
38     endif
39    
40 guez 110 call NF95_PUT_VAR(Ncid(Index), Varid = 1, values = Field, &
41 guez 109 start = (/1, Record(Index)/))
42    
43 guez 108 end subroutine WriteField1d
44    
45     !****************************************************************
46    
47     subroutine WriteField2d(name,Field)
48    
49 guez 113 use netcdf, only: nf90_sync
50    
51 guez 109 character(len=*), intent(in):: name
52     real, intent(in):: Field(:, :)
53 guez 108
54 guez 109 ! Local:
55 guez 113 integer index, status
56 guez 109
57     !-------------------------------------------
58    
59     Index=GetFieldIndex(name)
60    
61     if (Index==-1) then
62     call CreateNewField(name, shape(field))
63     Index=nbfield
64 guez 110 Record(Index) = 1
65 guez 109 else
66     Record(Index)=Record(Index)+1
67     endif
68    
69 guez 110 call NF95_PUT_VAR(Ncid(Index), Varid = 1, values = Field, &
70 guez 109 start = (/1, 1, Record(Index)/))
71 guez 113 status = nf90_sync(Ncid(Index))
72 guez 109
73 guez 108 end subroutine WriteField2d
74    
75     !****************************************************************
76    
77     subroutine WriteField3d(name,Field)
78    
79 guez 109 character(len=*), intent(in):: name
80     real, intent(in):: Field(:, :, :)
81 guez 108
82 guez 109 ! Local:
83     integer index
84    
85     !-------------------------------------------
86    
87     Index=GetFieldIndex(name)
88    
89     if (Index==-1) then
90     call CreateNewField(name, shape(field))
91     Index=nbfield
92 guez 110 Record(Index) = 1
93 guez 109 else
94     Record(Index)=Record(Index)+1
95     endif
96    
97 guez 110 call NF95_PUT_VAR(Ncid(Index), Varid = 1, values = Field, &
98 guez 109 start = (/1, 1, 1, Record(Index)/))
99    
100 guez 108 end subroutine WriteField3d
101    
102     end module WriteField_m

  ViewVC Help
Powered by ViewVC 1.1.21