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

Contents of /trunk/bibio/writefield.f

Parent Directory Parent Directory | Revision Log Revision Log


Revision 113 - (show annotations)
Thu Sep 18 19:56:46 2014 UTC (9 years, 7 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 module WriteField_m
2
3 use CreateNewField_m, only: CreateNewField, ncid, nbfield, MaxWriteField
4 use GetFieldIndex_m, only: GetFieldIndex
5 USE netcdf95, ONLY: nf95_put_var
6
7 implicit none
8
9 integer, save:: Record(MaxWriteField)
10
11 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 character(len=*), intent(in):: name
23 real, intent(in):: Field(:)
24
25 ! 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 Record(Index) = 1
36 else
37 Record(Index)=Record(Index)+1
38 endif
39
40 call NF95_PUT_VAR(Ncid(Index), Varid = 1, values = Field, &
41 start = (/1, Record(Index)/))
42
43 end subroutine WriteField1d
44
45 !****************************************************************
46
47 subroutine WriteField2d(name,Field)
48
49 use netcdf, only: nf90_sync
50
51 character(len=*), intent(in):: name
52 real, intent(in):: Field(:, :)
53
54 ! Local:
55 integer index, status
56
57 !-------------------------------------------
58
59 Index=GetFieldIndex(name)
60
61 if (Index==-1) then
62 call CreateNewField(name, shape(field))
63 Index=nbfield
64 Record(Index) = 1
65 else
66 Record(Index)=Record(Index)+1
67 endif
68
69 call NF95_PUT_VAR(Ncid(Index), Varid = 1, values = Field, &
70 start = (/1, 1, Record(Index)/))
71 status = nf90_sync(Ncid(Index))
72
73 end subroutine WriteField2d
74
75 !****************************************************************
76
77 subroutine WriteField3d(name,Field)
78
79 character(len=*), intent(in):: name
80 real, intent(in):: Field(:, :, :)
81
82 ! 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 Record(Index) = 1
93 else
94 Record(Index)=Record(Index)+1
95 endif
96
97 call NF95_PUT_VAR(Ncid(Index), Varid = 1, values = Field, &
98 start = (/1, 1, 1, Record(Index)/))
99
100 end subroutine WriteField3d
101
102 end module WriteField_m

  ViewVC Help
Powered by ViewVC 1.1.21