1 | ! |
---|
2 | ! $Id: write_field.f90 5270 2024-10-24 11:55:38Z fhourdin $ |
---|
3 | ! |
---|
4 | module write_field |
---|
5 | USE netcdf, ONLY: nf90_sync, nf90_put_var, nf90_enddef, nf90_def_dim, nf90_unlimited, & |
---|
6 | nf90_clobber, nf90_create, nf90_def_var, nf90_double |
---|
7 | USE strings_mod, ONLY: int2str |
---|
8 | IMPLICIT NONE; PRIVATE |
---|
9 | PUBLIC WriteField |
---|
10 | |
---|
11 | integer, parameter :: MaxWriteField = 100 |
---|
12 | integer, dimension(MaxWriteField),save :: FieldId |
---|
13 | integer, dimension(MaxWriteField),save :: FieldVarId |
---|
14 | integer, dimension(MaxWriteField),save :: FieldIndex |
---|
15 | character(len=255), dimension(MaxWriteField) :: FieldName |
---|
16 | |
---|
17 | integer,save :: NbField = 0 |
---|
18 | |
---|
19 | interface WriteField |
---|
20 | module procedure WriteField3d,WriteField2d,WriteField1d |
---|
21 | end interface WriteField |
---|
22 | contains |
---|
23 | |
---|
24 | function GetFieldIndex(name) |
---|
25 | implicit none |
---|
26 | integer :: GetFieldindex |
---|
27 | character(len=*) :: name |
---|
28 | |
---|
29 | character(len=255) :: TrueName |
---|
30 | integer :: i |
---|
31 | |
---|
32 | |
---|
33 | TrueName=TRIM(ADJUSTL(name)) |
---|
34 | |
---|
35 | GetFieldIndex=-1 |
---|
36 | do i=1,NbField |
---|
37 | if (TrueName==FieldName(i)) then |
---|
38 | GetFieldIndex=i |
---|
39 | exit |
---|
40 | endif |
---|
41 | enddo |
---|
42 | end function GetFieldIndex |
---|
43 | |
---|
44 | subroutine WriteField3d(name,Field) |
---|
45 | implicit none |
---|
46 | character(len=*) :: name |
---|
47 | real, dimension(:,:,:) :: Field |
---|
48 | integer, dimension(3) :: Dim |
---|
49 | |
---|
50 | Dim=shape(Field) |
---|
51 | call WriteField_gen(name,Field,Dim(1),Dim(2),Dim(3)) |
---|
52 | |
---|
53 | end subroutine WriteField3d |
---|
54 | |
---|
55 | subroutine WriteField2d(name,Field) |
---|
56 | implicit none |
---|
57 | character(len=*) :: name |
---|
58 | real, dimension(:,:) :: Field |
---|
59 | integer, dimension(2) :: Dim |
---|
60 | |
---|
61 | Dim=shape(Field) |
---|
62 | call WriteField_gen(name,Field,Dim(1),Dim(2),1) |
---|
63 | |
---|
64 | end subroutine WriteField2d |
---|
65 | |
---|
66 | subroutine WriteField1d(name,Field) |
---|
67 | implicit none |
---|
68 | character(len=*) :: name |
---|
69 | real, dimension(:) :: Field |
---|
70 | integer, dimension(1) :: Dim |
---|
71 | |
---|
72 | Dim=shape(Field) |
---|
73 | call WriteField_gen(name,Field,Dim(1),1,1) |
---|
74 | |
---|
75 | end subroutine WriteField1d |
---|
76 | |
---|
77 | subroutine WriteField_gen(name,Field,dimx,dimy,dimz) |
---|
78 | implicit none |
---|
79 | character(len=*) :: name |
---|
80 | integer :: dimx,dimy,dimz |
---|
81 | real,dimension(dimx,dimy,dimz) :: Field |
---|
82 | integer,dimension(dimx*dimy*dimz) :: ndex |
---|
83 | integer :: status |
---|
84 | integer :: index |
---|
85 | integer :: start(4) |
---|
86 | integer :: count(4) |
---|
87 | |
---|
88 | |
---|
89 | Index=GetFieldIndex(name) |
---|
90 | if (Index==-1) then |
---|
91 | call CreateNewField(name,dimx,dimy,dimz) |
---|
92 | Index=GetFieldIndex(name) |
---|
93 | else |
---|
94 | FieldIndex(Index)=FieldIndex(Index)+1. |
---|
95 | endif |
---|
96 | |
---|
97 | start(1)=1 |
---|
98 | start(2)=1 |
---|
99 | start(3)=1 |
---|
100 | start(4)=FieldIndex(Index) |
---|
101 | |
---|
102 | count(1)=dimx |
---|
103 | count(2)=dimy |
---|
104 | count(3)=dimz |
---|
105 | count(4)=1 |
---|
106 | |
---|
107 | status = nf90_put_var(FieldId(Index),FieldVarId(Index),Field,start,count) |
---|
108 | status = nf90_sync(FieldId(Index)) |
---|
109 | |
---|
110 | end subroutine WriteField_gen |
---|
111 | |
---|
112 | subroutine CreateNewField(name,dimx,dimy,dimz) |
---|
113 | implicit none |
---|
114 | character(len=*) :: name |
---|
115 | integer :: dimx,dimy,dimz |
---|
116 | integer :: TabDim(4) |
---|
117 | integer :: status |
---|
118 | |
---|
119 | |
---|
120 | NbField=NbField+1 |
---|
121 | FieldName(NbField)=TRIM(ADJUSTL(name)) |
---|
122 | FieldIndex(NbField)=1 |
---|
123 | |
---|
124 | |
---|
125 | status = nf90_create(TRIM(ADJUSTL(name))//'.nc', nf90_clobber, FieldId(NbField)) |
---|
126 | status = nf90_def_dim(FieldId(NbField),'X',dimx,TabDim(1)) |
---|
127 | status = nf90_def_dim(FieldId(NbField),'Y',dimy,TabDim(2)) |
---|
128 | status = nf90_def_dim(FieldId(NbField),'Z',dimz,TabDim(3)) |
---|
129 | status = nf90_def_dim(FieldId(NbField),'iter',nf90_unlimited,TabDim(4)) |
---|
130 | status = nf90_def_var(FieldId(NbField),FieldName(NbField),nf90_double,TabDim,FieldVarId(NbField)) |
---|
131 | status = nf90_enddef(FieldId(NbField)) |
---|
132 | |
---|
133 | end subroutine CreateNewField |
---|
134 | |
---|
135 | |
---|
136 | |
---|
137 | subroutine write_field1D(name,Field) |
---|
138 | implicit none |
---|
139 | |
---|
140 | integer, parameter :: MaxDim=1 |
---|
141 | character(len=*) :: name |
---|
142 | real, dimension(:) :: Field |
---|
143 | real, dimension(:),allocatable :: New_Field |
---|
144 | character(len=20) :: str |
---|
145 | integer, dimension(MaxDim) :: Dim |
---|
146 | integer :: i,nb |
---|
147 | integer, parameter :: id=10 |
---|
148 | integer, parameter :: NbCol=4 |
---|
149 | integer :: ColumnSize |
---|
150 | integer :: pos |
---|
151 | character(len=255) :: form |
---|
152 | character(len=255) :: MaxLen |
---|
153 | |
---|
154 | |
---|
155 | open(unit=id,file=name//'.field',form='formatted',status='replace') |
---|
156 | write (id,'("----- Field '//name//'",//)') |
---|
157 | Dim=shape(Field) |
---|
158 | MaxLen=int2str(len(trim(int2str(Dim(1))))) |
---|
159 | ColumnSize=20+6+3+len(trim(int2str(Dim(1)))) |
---|
160 | Nb=0 |
---|
161 | Pos=2 |
---|
162 | do i=1,Dim(1) |
---|
163 | nb=nb+1 |
---|
164 | |
---|
165 | if (MOD(nb,NbCol)==0) then |
---|
166 | form='(t'//trim(int2str(pos))// ',i'//trim(MaxLen) //'," ---> ",g22.16,/)' |
---|
167 | Pos=2 |
---|
168 | else |
---|
169 | form='(t'//trim(int2str(pos))// ',i'//trim(MaxLen) //'," ---> ",g22.16," | ",)' |
---|
170 | Pos=Pos+ColumnSize |
---|
171 | endif |
---|
172 | write (id,form,advance='no') i,Field(i) |
---|
173 | enddo |
---|
174 | |
---|
175 | close(id) |
---|
176 | |
---|
177 | end subroutine write_field1D |
---|
178 | |
---|
179 | subroutine write_field2D(name,Field) |
---|
180 | implicit none |
---|
181 | |
---|
182 | integer, parameter :: MaxDim=2 |
---|
183 | character(len=*) :: name |
---|
184 | real, dimension(:,:) :: Field |
---|
185 | real, dimension(:,:),allocatable :: New_Field |
---|
186 | character(len=20) :: str |
---|
187 | integer, dimension(MaxDim) :: Dim |
---|
188 | integer :: i,j,nb |
---|
189 | integer, parameter :: id=10 |
---|
190 | integer, parameter :: NbCol=4 |
---|
191 | integer :: ColumnSize |
---|
192 | integer :: pos,offset |
---|
193 | character(len=255) :: form |
---|
194 | character(len=255) :: spacing |
---|
195 | |
---|
196 | open(unit=id,file=name//'.field',form='formatted',status='replace') |
---|
197 | write (id,'("----- Field '//name//'",//)') |
---|
198 | |
---|
199 | Dim=shape(Field) |
---|
200 | offset=len(trim(int2str(Dim(1))))+len(trim(int2str(Dim(2))))+3 |
---|
201 | ColumnSize=20+6+3+offset |
---|
202 | |
---|
203 | spacing='(t2,"'//repeat('-',ColumnSize*NbCol)//'")' |
---|
204 | |
---|
205 | do i=1,Dim(2) |
---|
206 | nb=0 |
---|
207 | Pos=2 |
---|
208 | do j=1,Dim(1) |
---|
209 | nb=nb+1 |
---|
210 | |
---|
211 | if (MOD(nb,NbCol)==0) then |
---|
212 | form='(t'//trim(int2str(pos))// & |
---|
213 | ',"('//trim(int2str(j))//',' & |
---|
214 | //trim(int2str(i))//')",t' & |
---|
215 | //trim(int2str(pos+offset)) & |
---|
216 | //'," ---> ",g22.16,/)' |
---|
217 | Pos=2 |
---|
218 | else |
---|
219 | form='(t'//trim(int2str(pos))// & |
---|
220 | ',"('//trim(int2str(j))//',' & |
---|
221 | //trim(int2str(i))//')",t' & |
---|
222 | //trim(int2str(pos+offset)) & |
---|
223 | //'," ---> ",g22.16," | ")' |
---|
224 | Pos=Pos+ColumnSize |
---|
225 | endif |
---|
226 | write (id,form,advance='no') Field(j,i) |
---|
227 | enddo |
---|
228 | if (MOD(nb,NbCol)==0) then |
---|
229 | write (id,spacing) |
---|
230 | else |
---|
231 | write (id,'("")') |
---|
232 | write (id,spacing) |
---|
233 | endif |
---|
234 | enddo |
---|
235 | |
---|
236 | end subroutine write_field2D |
---|
237 | |
---|
238 | subroutine write_field3D(name,Field) |
---|
239 | implicit none |
---|
240 | |
---|
241 | integer, parameter :: MaxDim=3 |
---|
242 | character(len=*) :: name |
---|
243 | real, dimension(:,:,:) :: Field |
---|
244 | real, dimension(:,:,:),allocatable :: New_Field |
---|
245 | integer, dimension(MaxDim) :: Dim |
---|
246 | integer :: i,j,k,nb |
---|
247 | integer, parameter :: id=10 |
---|
248 | integer, parameter :: NbCol=4 |
---|
249 | integer :: ColumnSize |
---|
250 | integer :: pos,offset |
---|
251 | character(len=255) :: form |
---|
252 | character(len=255) :: spacing |
---|
253 | |
---|
254 | open(unit=id,file=name//'.field',form='formatted',status='replace') |
---|
255 | write (id,'("----- Field '//name//'"//)') |
---|
256 | |
---|
257 | Dim=shape(Field) |
---|
258 | offset=len(trim(int2str(Dim(1))))+len(trim(int2str(Dim(2))))+len(trim(int2str(Dim(3))))+4 |
---|
259 | ColumnSize=22+6+3+offset |
---|
260 | |
---|
261 | ! open(unit=id,file=name,form=formatted |
---|
262 | |
---|
263 | spacing='(t2,"'//repeat('-',ColumnSize*NbCol)//'")' |
---|
264 | |
---|
265 | do i=1,Dim(3) |
---|
266 | |
---|
267 | do j=1,Dim(2) |
---|
268 | nb=0 |
---|
269 | Pos=2 |
---|
270 | |
---|
271 | do k=1,Dim(1) |
---|
272 | nb=nb+1 |
---|
273 | |
---|
274 | if (MOD(nb,NbCol)==0) then |
---|
275 | form='(t'//trim(int2str(pos))// & |
---|
276 | ',"('//trim(int2str(k))//',' & |
---|
277 | //trim(int2str(j))//',' & |
---|
278 | //trim(int2str(i))//')",t' & |
---|
279 | //trim(int2str(pos+offset)) & |
---|
280 | //'," ---> ",g22.16,/)' |
---|
281 | Pos=2 |
---|
282 | else |
---|
283 | form='(t'//trim(int2str(pos))// & |
---|
284 | ',"('//trim(int2str(k))//',' & |
---|
285 | //trim(int2str(j))//',' & |
---|
286 | //trim(int2str(i))//')",t' & |
---|
287 | //trim(int2str(pos+offset)) & |
---|
288 | //'," ---> ",g22.16," | ")' |
---|
289 | ! d�pent de l'impl�mention, sur compaq, c'est necessaire |
---|
290 | ! Pos=Pos+ColumnSize |
---|
291 | endif |
---|
292 | write (id,form,advance='no') Field(k,j,i) |
---|
293 | enddo |
---|
294 | if (MOD(nb,NbCol)==0) then |
---|
295 | write (id,spacing) |
---|
296 | else |
---|
297 | write (id,'("")') |
---|
298 | write (id,spacing) |
---|
299 | endif |
---|
300 | enddo |
---|
301 | write (id,spacing) |
---|
302 | enddo |
---|
303 | |
---|
304 | close(id) |
---|
305 | |
---|
306 | end subroutine write_field3D |
---|
307 | |
---|
308 | end module write_field |
---|
309 | |
---|