1 | |
---|
2 | ! $Id$ |
---|
3 | |
---|
4 | module Bands |
---|
5 | USE parallel_lmdz |
---|
6 | integer, parameter :: bands_caldyn=1 |
---|
7 | integer, parameter :: bands_vanleer=2 |
---|
8 | integer, parameter :: bands_dissip=3 |
---|
9 | |
---|
10 | INTEGER,dimension(:),allocatable :: jj_Nb_Caldyn |
---|
11 | INTEGER,dimension(:),allocatable :: jj_Nb_vanleer |
---|
12 | INTEGER,dimension(:),allocatable :: jj_Nb_vanleer2 |
---|
13 | INTEGER,dimension(:),allocatable :: jj_Nb_dissip |
---|
14 | INTEGER,dimension(:),allocatable :: jj_Nb_physic |
---|
15 | INTEGER,dimension(:),allocatable :: jj_Nb_physic_bis |
---|
16 | |
---|
17 | TYPE(distrib),SAVE,TARGET :: distrib_Caldyn |
---|
18 | TYPE(distrib),SAVE,TARGET :: distrib_vanleer |
---|
19 | TYPE(distrib),SAVE,TARGET :: distrib_vanleer2 |
---|
20 | TYPE(distrib),SAVE,TARGET :: distrib_dissip |
---|
21 | TYPE(distrib),SAVE,TARGET :: distrib_physic |
---|
22 | TYPE(distrib),SAVE,TARGET :: distrib_physic_bis |
---|
23 | |
---|
24 | INTEGER,dimension(:),allocatable :: distrib_phys |
---|
25 | |
---|
26 | contains |
---|
27 | |
---|
28 | SUBROUTINE AllocateBands |
---|
29 | USE parallel_lmdz |
---|
30 | IMPLICIT NONE |
---|
31 | |
---|
32 | allocate(jj_Nb_Caldyn(0:MPI_Size-1)) |
---|
33 | allocate(jj_Nb_vanleer(0:MPI_Size-1)) |
---|
34 | allocate(jj_Nb_vanleer2(0:MPI_Size-1)) |
---|
35 | allocate(jj_Nb_dissip(0:MPI_Size-1)) |
---|
36 | allocate(jj_Nb_physic(0:MPI_Size-1)) |
---|
37 | allocate(jj_Nb_physic_bis(0:MPI_Size-1)) |
---|
38 | allocate(distrib_phys(0:MPI_Size-1)) |
---|
39 | |
---|
40 | END SUBROUTINE AllocateBands |
---|
41 | |
---|
42 | SUBROUTINE Read_distrib |
---|
43 | USE parallel_lmdz |
---|
44 | IMPLICIT NONE |
---|
45 | |
---|
46 | include "dimensions.h" |
---|
47 | INTEGER :: i,j |
---|
48 | character (len=4) :: siim,sjjm,sllm,sproc |
---|
49 | character (len=255) :: filename |
---|
50 | INTEGER :: unit_number=10 |
---|
51 | INTEGER :: ierr |
---|
52 | |
---|
53 | CALL AllocateBands |
---|
54 | WRITE(siim,'(i3)') iim |
---|
55 | WRITE(sjjm,'(i3)') jjm |
---|
56 | WRITE(sllm,'(i3)') llm |
---|
57 | WRITE(sproc,'(i3)') mpi_size |
---|
58 | filename='Bands_'//TRIM(ADJUSTL(siim))//'x'//TRIM(ADJUSTL(sjjm))//'x'//TRIM(ADJUSTL(sllm))//'_' & |
---|
59 | //TRIM(ADJUSTL(sproc))//'prc.dat' |
---|
60 | |
---|
61 | OPEN(UNIT=unit_number,FILE=trim(filename),STATUS='old',FORM='formatted',IOSTAT=ierr) |
---|
62 | |
---|
63 | if (ierr==0) THEN |
---|
64 | do i=0,mpi_size-1 |
---|
65 | read (unit_number,*) j,jj_nb_caldyn(i) |
---|
66 | enddo |
---|
67 | |
---|
68 | do i=0,mpi_size-1 |
---|
69 | read (unit_number,*) j,jj_nb_vanleer(i) |
---|
70 | enddo |
---|
71 | |
---|
72 | do i=0,mpi_size-1 |
---|
73 | read (unit_number,*) j,jj_nb_dissip(i) |
---|
74 | enddo |
---|
75 | |
---|
76 | do i=0,mpi_size-1 |
---|
77 | read (unit_number,*) j,distrib_phys(i) |
---|
78 | enddo |
---|
79 | |
---|
80 | CLOSE(unit_number) |
---|
81 | |
---|
82 | else |
---|
83 | do i=0,mpi_size-1 |
---|
84 | jj_nb_caldyn(i)=(jjm+1)/mpi_size |
---|
85 | if (i<MOD(jjm+1,mpi_size)) jj_nb_caldyn(i)=jj_nb_caldyn(i)+1 |
---|
86 | enddo |
---|
87 | |
---|
88 | jj_nb_vanleer(:)=jj_nb_caldyn(:) |
---|
89 | jj_nb_dissip(:)=jj_nb_caldyn(:) |
---|
90 | |
---|
91 | do i=0,mpi_size-1 |
---|
92 | distrib_phys(i)=(iim*(jjm-1)+2)/mpi_size |
---|
93 | IF (i<MOD(iim*(jjm-1)+2,mpi_size)) distrib_phys(i)=distrib_phys(i)+1 |
---|
94 | enddo |
---|
95 | endif |
---|
96 | |
---|
97 | ! distrib_phys(:)=jj_nb_caldyn(:)*iim |
---|
98 | ! distrib_phys(0) = distrib_phys(0) - (iim-1) |
---|
99 | ! distrib_phys(mpi_size-1) = distrib_phys(mpi_size-1) - (iim-1) |
---|
100 | |
---|
101 | END SUBROUTINE Read_distrib |
---|
102 | |
---|
103 | |
---|
104 | SUBROUTINE Set_Bands |
---|
105 | USE parallel_lmdz |
---|
106 | IMPLICIT NONE |
---|
107 | INCLUDE 'dimensions.h' |
---|
108 | INTEGER :: i, ij |
---|
109 | INTEGER :: jj_para_begin(0:mpi_size-1) |
---|
110 | INTEGER :: jj_para_end(0:mpi_size-1) |
---|
111 | |
---|
112 | do i=0,mpi_size-1 |
---|
113 | jj_nb_vanleer2(i)=(jjm+1)/mpi_size |
---|
114 | if (i<MOD(jjm+1,mpi_size)) jj_nb_vanleer2(i)=jj_nb_vanleer2(i)+1 |
---|
115 | enddo |
---|
116 | |
---|
117 | jj_para_begin(0)=1 |
---|
118 | ij=distrib_phys(0)+iim-1 |
---|
119 | jj_para_end(0)=((ij-1)/iim)+1 |
---|
120 | |
---|
121 | DO i=1,mpi_Size-1 |
---|
122 | ij=ij+1 |
---|
123 | jj_para_begin(i)=((ij-1)/iim)+1 |
---|
124 | ij=ij+distrib_phys(i)-1 |
---|
125 | jj_para_end(i)=((ij-1)/iim)+1 |
---|
126 | ENDDO |
---|
127 | |
---|
128 | do i=0,MPI_Size-1 |
---|
129 | jj_Nb_physic(i)=jj_para_end(i)-jj_para_begin(i)+1 |
---|
130 | if (i/=0) THEN |
---|
131 | if (jj_para_begin(i)==jj_para_end(i-1)) THEN |
---|
132 | jj_Nb_physic(i-1)=jj_Nb_physic(i-1)-1 |
---|
133 | endif |
---|
134 | endif |
---|
135 | enddo |
---|
136 | |
---|
137 | do i=0,MPI_Size-1 |
---|
138 | jj_Nb_physic_bis(i)=jj_para_end(i)-jj_para_begin(i)+1 |
---|
139 | if (i/=0) THEN |
---|
140 | if (jj_para_begin(i)==jj_para_end(i-1)) THEN |
---|
141 | jj_Nb_physic_bis(i)=jj_Nb_physic_bis(i)-1 |
---|
142 | else |
---|
143 | jj_Nb_physic_bis(i-1)=jj_Nb_physic_bis(i-1)+1 |
---|
144 | jj_Nb_physic_bis(i)=jj_Nb_physic_bis(i)-1 |
---|
145 | endif |
---|
146 | endif |
---|
147 | enddo |
---|
148 | |
---|
149 | CALL create_distrib(jj_Nb_Caldyn,distrib_caldyn) |
---|
150 | CALL create_distrib(jj_Nb_vanleer,distrib_vanleer) |
---|
151 | CALL create_distrib(jj_Nb_vanleer2,distrib_vanleer2) |
---|
152 | CALL create_distrib(jj_Nb_dissip,distrib_dissip) |
---|
153 | CALL create_distrib(jj_Nb_physic,distrib_physic) |
---|
154 | CALL create_distrib(jj_Nb_physic_bis,distrib_physic_bis) |
---|
155 | |
---|
156 | distrib_physic_bis%jjb_u=distrib_physic%jjb_u |
---|
157 | distrib_physic_bis%jje_u=distrib_physic%jje_u |
---|
158 | distrib_physic_bis%jjnb_u=distrib_physic%jjnb_u |
---|
159 | |
---|
160 | distrib_physic_bis%ijb_u=distrib_physic%ijb_u |
---|
161 | distrib_physic_bis%ije_u=distrib_physic%ije_u |
---|
162 | distrib_physic_bis%ijnb_u=distrib_physic%ijnb_u |
---|
163 | |
---|
164 | distrib_physic_bis%jjb_v=distrib_physic%jjb_v |
---|
165 | distrib_physic_bis%jje_v=distrib_physic%jje_v |
---|
166 | distrib_physic_bis%jjnb_v=distrib_physic%jjnb_v |
---|
167 | |
---|
168 | distrib_physic_bis%ijb_v=distrib_physic%ijb_v |
---|
169 | distrib_physic_bis%ije_v=distrib_physic%ije_v |
---|
170 | distrib_physic_bis%ijnb_v=distrib_physic%ijnb_v |
---|
171 | |
---|
172 | END SUBROUTINE Set_Bands |
---|
173 | |
---|
174 | |
---|
175 | SUBROUTINE AdjustBands_caldyn(new_dist) |
---|
176 | use times |
---|
177 | USE parallel_lmdz |
---|
178 | IMPLICIT NONE |
---|
179 | TYPE(distrib),INTENT(INOUT) :: new_dist |
---|
180 | |
---|
181 | REAL :: minvalue,maxvalue |
---|
182 | INTEGER :: min_proc,max_proc |
---|
183 | INTEGER :: i,j |
---|
184 | real,allocatable,dimension(:) :: value |
---|
185 | integer,allocatable,dimension(:) :: index |
---|
186 | REAL :: tmpvalue |
---|
187 | INTEGER :: tmpindex |
---|
188 | |
---|
189 | allocate(value(0:mpi_size-1)) |
---|
190 | allocate(index(0:mpi_size-1)) |
---|
191 | |
---|
192 | |
---|
193 | CALL allgather_timer_average |
---|
194 | |
---|
195 | do i=0,mpi_size-1 |
---|
196 | value(i)=timer_average(jj_nb_caldyn(i),timer_caldyn,i) |
---|
197 | index(i)=i |
---|
198 | enddo |
---|
199 | |
---|
200 | do i=0,mpi_size-2 |
---|
201 | do j=i+1,mpi_size-1 |
---|
202 | if (value(i)>value(j)) THEN |
---|
203 | tmpvalue=value(i) |
---|
204 | value(i)=value(j) |
---|
205 | value(j)=tmpvalue |
---|
206 | |
---|
207 | tmpindex=index(i) |
---|
208 | index(i)=index(j) |
---|
209 | index(j)=tmpindex |
---|
210 | endif |
---|
211 | enddo |
---|
212 | enddo |
---|
213 | |
---|
214 | maxvalue=value(mpi_size-1) |
---|
215 | max_proc=index(mpi_size-1) |
---|
216 | |
---|
217 | do i=0,mpi_size-2 |
---|
218 | minvalue=value(i) |
---|
219 | min_proc=index(i) |
---|
220 | if (jj_nb_caldyn(max_proc)>2) THEN |
---|
221 | if (timer_iteration(jj_nb_caldyn(min_proc)+1,timer_caldyn,min_proc)<=1 ) THEN |
---|
222 | jj_nb_caldyn(min_proc)=jj_nb_caldyn(min_proc)+1 |
---|
223 | jj_nb_caldyn(max_proc)=jj_nb_caldyn(max_proc)-1 |
---|
224 | exit |
---|
225 | else |
---|
226 | if (timer_average(jj_nb_caldyn(min_proc)+1,timer_caldyn,min_proc) & |
---|
227 | -timer_delta(jj_nb_caldyn(min_proc)+1,timer_caldyn,min_proc) < maxvalue) THEN |
---|
228 | jj_nb_caldyn(min_proc)=jj_nb_caldyn(min_proc)+1 |
---|
229 | jj_nb_caldyn(max_proc)=jj_nb_caldyn(max_proc)-1 |
---|
230 | exit |
---|
231 | endif |
---|
232 | endif |
---|
233 | endif |
---|
234 | enddo |
---|
235 | |
---|
236 | deallocate(value) |
---|
237 | deallocate(index) |
---|
238 | CALL create_distrib(jj_nb_caldyn,new_dist) |
---|
239 | |
---|
240 | END SUBROUTINE AdjustBands_caldyn |
---|
241 | |
---|
242 | SUBROUTINE AdjustBands_vanleer(new_dist) |
---|
243 | use times |
---|
244 | USE parallel_lmdz |
---|
245 | IMPLICIT NONE |
---|
246 | TYPE(distrib),INTENT(INOUT) :: new_dist |
---|
247 | |
---|
248 | REAL :: minvalue,maxvalue |
---|
249 | INTEGER :: min_proc,max_proc |
---|
250 | INTEGER :: i,j |
---|
251 | real,allocatable,dimension(:) :: value |
---|
252 | integer,allocatable,dimension(:) :: index |
---|
253 | REAL :: tmpvalue |
---|
254 | INTEGER :: tmpindex |
---|
255 | |
---|
256 | allocate(value(0:mpi_size-1)) |
---|
257 | allocate(index(0:mpi_size-1)) |
---|
258 | |
---|
259 | |
---|
260 | CALL allgather_timer_average |
---|
261 | |
---|
262 | do i=0,mpi_size-1 |
---|
263 | value(i)=timer_average(jj_nb_vanleer(i),timer_vanleer,i) |
---|
264 | index(i)=i |
---|
265 | enddo |
---|
266 | |
---|
267 | do i=0,mpi_size-2 |
---|
268 | do j=i+1,mpi_size-1 |
---|
269 | if (value(i)>value(j)) THEN |
---|
270 | tmpvalue=value(i) |
---|
271 | value(i)=value(j) |
---|
272 | value(j)=tmpvalue |
---|
273 | |
---|
274 | tmpindex=index(i) |
---|
275 | index(i)=index(j) |
---|
276 | index(j)=tmpindex |
---|
277 | endif |
---|
278 | enddo |
---|
279 | enddo |
---|
280 | |
---|
281 | maxvalue=value(mpi_size-1) |
---|
282 | max_proc=index(mpi_size-1) |
---|
283 | |
---|
284 | do i=0,mpi_size-2 |
---|
285 | minvalue=value(i) |
---|
286 | min_proc=index(i) |
---|
287 | |
---|
288 | if (jj_nb_vanleer(max_proc)>2) THEN |
---|
289 | if (timer_average(jj_nb_vanleer(min_proc)+1,timer_vanleer,min_proc)==0. .or. & |
---|
290 | timer_average(jj_nb_vanleer(max_proc)-1,timer_vanleer,max_proc)==0.) THEN |
---|
291 | jj_nb_vanleer(min_proc)=jj_nb_vanleer(min_proc)+1 |
---|
292 | jj_nb_vanleer(max_proc)=jj_nb_vanleer(max_proc)-1 |
---|
293 | exit |
---|
294 | else |
---|
295 | if (timer_average(jj_nb_vanleer(min_proc)+1,timer_vanleer,min_proc) < maxvalue) THEN |
---|
296 | jj_nb_vanleer(min_proc)=jj_nb_vanleer(min_proc)+1 |
---|
297 | jj_nb_vanleer(max_proc)=jj_nb_vanleer(max_proc)-1 |
---|
298 | exit |
---|
299 | endif |
---|
300 | endif |
---|
301 | endif |
---|
302 | enddo |
---|
303 | |
---|
304 | deallocate(value) |
---|
305 | deallocate(index) |
---|
306 | |
---|
307 | CALL create_distrib(jj_nb_vanleer,new_dist) |
---|
308 | |
---|
309 | END SUBROUTINE AdjustBands_vanleer |
---|
310 | |
---|
311 | SUBROUTINE AdjustBands_dissip(new_dist) |
---|
312 | use times |
---|
313 | USE parallel_lmdz |
---|
314 | IMPLICIT NONE |
---|
315 | TYPE(distrib),INTENT(INOUT) :: new_dist |
---|
316 | |
---|
317 | REAL :: minvalue,maxvalue |
---|
318 | INTEGER :: min_proc,max_proc |
---|
319 | INTEGER :: i,j |
---|
320 | real,allocatable,dimension(:) :: value |
---|
321 | integer,allocatable,dimension(:) :: index |
---|
322 | REAL :: tmpvalue |
---|
323 | INTEGER :: tmpindex |
---|
324 | |
---|
325 | allocate(value(0:mpi_size-1)) |
---|
326 | allocate(index(0:mpi_size-1)) |
---|
327 | |
---|
328 | |
---|
329 | CALL allgather_timer_average |
---|
330 | |
---|
331 | do i=0,mpi_size-1 |
---|
332 | value(i)=timer_average(jj_nb_dissip(i),timer_dissip,i) |
---|
333 | index(i)=i |
---|
334 | enddo |
---|
335 | |
---|
336 | do i=0,mpi_size-2 |
---|
337 | do j=i+1,mpi_size-1 |
---|
338 | if (value(i)>value(j)) THEN |
---|
339 | tmpvalue=value(i) |
---|
340 | value(i)=value(j) |
---|
341 | value(j)=tmpvalue |
---|
342 | |
---|
343 | tmpindex=index(i) |
---|
344 | index(i)=index(j) |
---|
345 | index(j)=tmpindex |
---|
346 | endif |
---|
347 | enddo |
---|
348 | enddo |
---|
349 | |
---|
350 | maxvalue=value(mpi_size-1) |
---|
351 | max_proc=index(mpi_size-1) |
---|
352 | |
---|
353 | do i=0,mpi_size-2 |
---|
354 | minvalue=value(i) |
---|
355 | min_proc=index(i) |
---|
356 | |
---|
357 | if (jj_nb_dissip(max_proc)>3) THEN |
---|
358 | if (timer_iteration(jj_nb_dissip(min_proc)+1,timer_dissip,min_proc)<=1) THEN |
---|
359 | jj_nb_dissip(min_proc)=jj_nb_dissip(min_proc)+1 |
---|
360 | jj_nb_dissip(max_proc)=jj_nb_dissip(max_proc)-1 |
---|
361 | exit |
---|
362 | else |
---|
363 | if (timer_average(jj_nb_dissip(min_proc)+1,timer_dissip,min_proc) & |
---|
364 | - timer_delta(jj_nb_dissip(min_proc)+1,timer_dissip,min_proc) < maxvalue) THEN |
---|
365 | jj_nb_dissip(min_proc)=jj_nb_dissip(min_proc)+1 |
---|
366 | jj_nb_dissip(max_proc)=jj_nb_dissip(max_proc)-1 |
---|
367 | exit |
---|
368 | endif |
---|
369 | endif |
---|
370 | endif |
---|
371 | enddo |
---|
372 | |
---|
373 | deallocate(value) |
---|
374 | deallocate(index) |
---|
375 | |
---|
376 | CALL create_distrib(jj_nb_dissip,new_dist) |
---|
377 | |
---|
378 | END SUBROUTINE AdjustBands_dissip |
---|
379 | |
---|
380 | SUBROUTINE AdjustBands_physic |
---|
381 | use times |
---|
382 | |
---|
383 | ! Ehouarn: what follows is only related to // physics |
---|
384 | USE lmdz_cppkeys_wrapper, ONLY: CPPKEY_PHYS |
---|
385 | USE lmdz_phys_para, ONLY: klon_mpi_para_nb |
---|
386 | |
---|
387 | USE parallel_lmdz |
---|
388 | IMPLICIT NONE |
---|
389 | |
---|
390 | INTEGER :: i,Index |
---|
391 | real,allocatable,dimension(:) :: value |
---|
392 | integer,allocatable,dimension(:) :: Inc |
---|
393 | REAL :: medium |
---|
394 | INTEGER :: NbTot,sgn |
---|
395 | |
---|
396 | allocate(value(0:mpi_size-1)) |
---|
397 | allocate(Inc(0:mpi_size-1)) |
---|
398 | |
---|
399 | |
---|
400 | CALL allgather_timer_average |
---|
401 | |
---|
402 | medium=0 |
---|
403 | do i=0,mpi_size-1 |
---|
404 | value(i)=timer_average(jj_nb_physic(i),timer_physic,i) |
---|
405 | medium=medium+value(i) |
---|
406 | enddo |
---|
407 | |
---|
408 | medium=medium/mpi_size |
---|
409 | NbTot=0 |
---|
410 | IF (CPPKEY_PHYS) THEN |
---|
411 | do i=0,mpi_size-1 |
---|
412 | Inc(i)=nint(klon_mpi_para_nb(i)*(medium-value(i))/value(i)) |
---|
413 | NbTot=NbTot+Inc(i) |
---|
414 | enddo |
---|
415 | |
---|
416 | if (NbTot>=0) THEN |
---|
417 | Sgn=1 |
---|
418 | else |
---|
419 | Sgn=-1 |
---|
420 | NbTot=-NbTot |
---|
421 | endif |
---|
422 | |
---|
423 | Index=0 |
---|
424 | do i=1,NbTot |
---|
425 | Inc(Index)=Inc(Index)-Sgn |
---|
426 | Index=Index+1 |
---|
427 | if (Index>mpi_size-1) Index=0 |
---|
428 | enddo |
---|
429 | |
---|
430 | do i=0,mpi_size-1 |
---|
431 | distrib_phys(i)=klon_mpi_para_nb(i)+inc(i) |
---|
432 | enddo |
---|
433 | END IF |
---|
434 | |
---|
435 | END SUBROUTINE AdjustBands_physic |
---|
436 | |
---|
437 | SUBROUTINE WriteBands |
---|
438 | USE parallel_lmdz |
---|
439 | IMPLICIT NONE |
---|
440 | include "dimensions.h" |
---|
441 | |
---|
442 | INTEGER :: i,j |
---|
443 | character (len=4) :: siim,sjjm,sllm,sproc |
---|
444 | character (len=255) :: filename |
---|
445 | INTEGER :: unit_number=10 |
---|
446 | INTEGER :: ierr |
---|
447 | |
---|
448 | WRITE(siim,'(i3)') iim |
---|
449 | WRITE(sjjm,'(i3)') jjm |
---|
450 | WRITE(sllm,'(i3)') llm |
---|
451 | WRITE(sproc,'(i3)') mpi_size |
---|
452 | |
---|
453 | filename='Bands_'//TRIM(ADJUSTL(siim))//'x'//TRIM(ADJUSTL(sjjm))//'x'//TRIM(ADJUSTL(sllm))//'_' & |
---|
454 | //TRIM(ADJUSTL(sproc))//'prc.dat' |
---|
455 | |
---|
456 | OPEN(UNIT=unit_number,FILE=trim(filename),STATUS='replace',FORM='formatted',IOSTAT=ierr) |
---|
457 | |
---|
458 | if (ierr==0) THEN |
---|
459 | ! write (unit_number,*) '*** Bandes caldyn ***' |
---|
460 | do i=0,mpi_size-1 |
---|
461 | write (unit_number,*) i,jj_nb_caldyn(i) |
---|
462 | enddo |
---|
463 | |
---|
464 | ! write (unit_number,*) '*** Bandes vanleer ***' |
---|
465 | do i=0,mpi_size-1 |
---|
466 | write (unit_number,*) i,jj_nb_vanleer(i) |
---|
467 | enddo |
---|
468 | |
---|
469 | ! write (unit_number,*) '*** Bandes dissip ***' |
---|
470 | do i=0,mpi_size-1 |
---|
471 | write (unit_number,*) i,jj_nb_dissip(i) |
---|
472 | enddo |
---|
473 | |
---|
474 | do i=0,mpi_size-1 |
---|
475 | write (unit_number,*) i,distrib_phys(i) |
---|
476 | enddo |
---|
477 | |
---|
478 | CLOSE(unit_number) |
---|
479 | else |
---|
480 | print *,'probleme lors de l ecriture des bandes' |
---|
481 | endif |
---|
482 | |
---|
483 | END SUBROUTINE WriteBands |
---|
484 | |
---|
485 | end module Bands |
---|
486 | |
---|
487 | |
---|
488 | |
---|