LMDZ
mod_phys_lmdz_omp_transfert.F90
Go to the documentation of this file.
1 !
2 !$Id$
3 !
5 
6  PRIVATE
7 
8  INTEGER,PARAMETER :: grow_factor=1.5
9  INTEGER,PARAMETER :: size_min=1024
10 
11  CHARACTER(LEN=size_min),SAVE :: buffer_c
12 ! INTEGER,SAVE :: size_c=0
13  INTEGER,SAVE,ALLOCATABLE,DIMENSION(:) :: buffer_i
14  INTEGER,SAVE :: size_i=0
15  REAL,SAVE,ALLOCATABLE,DIMENSION(:) :: buffer_r
16  INTEGER,SAVE :: size_r=0
17  LOGICAL,SAVE,ALLOCATABLE,DIMENSION(:) :: buffer_l
18  INTEGER,SAVE :: size_l=0
19 
20 
21 
22 
23  INTERFACE bcast_omp
24  MODULE PROCEDURE bcast_omp_c, &
28  END INTERFACE
29 
30  INTERFACE scatter_omp
34  END INTERFACE
35 
36 
37  INTERFACE gather_omp
41  END INTERFACE
42 
43 
44  INTERFACE reduce_sum_omp
47  END INTERFACE
48 
49 
51 
52 CONTAINS
53 
54  SUBROUTINE omp_barrier
55  IMPLICIT NONE
56 
57 !$OMP BARRIER
58 
59  END SUBROUTINE omp_barrier
60 
61  SUBROUTINE check_buffer_i(buff_size)
62  IMPLICIT NONE
63  INTEGER :: buff_size
64 
65 !$OMP BARRIER
66 !$OMP MASTER
67  IF (buff_size>size_i) THEN
68  IF (ALLOCATED(buffer_i)) DEALLOCATE(buffer_i)
69  size_i=max(size_min,int(grow_factor*buff_size))
70  ALLOCATE(buffer_i(size_i))
71  ENDIF
72 !$OMP END MASTER
73 !$OMP BARRIER
74 
75  END SUBROUTINE check_buffer_i
76 
77  SUBROUTINE check_buffer_r(buff_size)
78  IMPLICIT NONE
79  INTEGER :: buff_size
80 
81 !$OMP BARRIER
82 !$OMP MASTER
83  IF (buff_size>size_r) THEN
84  IF (ALLOCATED(buffer_r)) DEALLOCATE(buffer_r)
85  size_r=max(size_min,int(grow_factor*buff_size))
86  ALLOCATE(buffer_r(size_r))
87  ENDIF
88 !$OMP END MASTER
89 !$OMP BARRIER
90 
91  END SUBROUTINE check_buffer_r
92 
93  SUBROUTINE check_buffer_l(buff_size)
94  IMPLICIT NONE
95  INTEGER :: buff_size
96 
97 !$OMP BARRIER
98 !$OMP MASTER
99  IF (buff_size>size_l) THEN
100  IF (ALLOCATED(buffer_l)) DEALLOCATE(buffer_l)
101  size_l=max(size_min,int(grow_factor*buff_size))
102  ALLOCATE(buffer_l(size_l))
103  ENDIF
104 !$OMP END MASTER
105 !$OMP BARRIER
106 
107  END SUBROUTINE check_buffer_l
108 
109 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
110 !! Definition des Broadcast --> 4D !!
111 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
112 
113 !! -- Les chaine de charactère -- !!
114 
115  SUBROUTINE bcast_omp_c(var)
116  IMPLICIT NONE
117  CHARACTER(LEN=*),INTENT(INOUT) :: Var
118 
119  CALL bcast_omp_cgen(var,len(var),buffer_c)
120 
121  END SUBROUTINE bcast_omp_c
122 
123 !! -- Les entiers -- !!
124 
125  SUBROUTINE bcast_omp_i(var)
126  IMPLICIT NONE
127  INTEGER,INTENT(INOUT) :: Var
128  INTEGER :: Var_tmp(1)
129 
130  var_tmp(1)=var
131  CALL check_buffer_i(1)
132  CALL bcast_omp_igen(var_tmp,1,buffer_i)
133  var=var_tmp(1)
134 
135  END SUBROUTINE bcast_omp_i
136 
137 
138  SUBROUTINE bcast_omp_i1(var)
139  IMPLICIT NONE
140  INTEGER,INTENT(INOUT) :: Var(:)
141 
142  CALL check_buffer_i(size(var))
143  CALL bcast_omp_igen(var,size(var),buffer_i)
144 
145  END SUBROUTINE bcast_omp_i1
146 
147 
148  SUBROUTINE bcast_omp_i2(var)
149  IMPLICIT NONE
150  INTEGER,INTENT(INOUT) :: Var(:,:)
151 
152  CALL check_buffer_i(size(var))
153  CALL bcast_omp_igen(var,size(var),buffer_i)
154 
155  END SUBROUTINE bcast_omp_i2
156 
157 
158  SUBROUTINE bcast_omp_i3(var)
159  IMPLICIT NONE
160  INTEGER,INTENT(INOUT) :: Var(:,:,:)
161 
162  CALL check_buffer_i(size(var))
163  CALL bcast_omp_igen(var,size(var),buffer_i)
164 
165  END SUBROUTINE bcast_omp_i3
166 
167 
168  SUBROUTINE bcast_omp_i4(var)
169  IMPLICIT NONE
170  INTEGER,INTENT(INOUT) :: Var(:,:,:,:)
171 
172  CALL check_buffer_i(size(var))
173  CALL bcast_omp_igen(var,size(var),buffer_i)
174 
175  END SUBROUTINE bcast_omp_i4
176 
177 
178 !! -- Les reels -- !!
179 
180  SUBROUTINE bcast_omp_r(var)
181  IMPLICIT NONE
182  REAL,INTENT(INOUT) :: Var
183  REAL :: Var_tmp(1)
184 
185  var_tmp(1)=var
186  CALL check_buffer_r(1)
187  CALL bcast_omp_rgen(var_tmp,1,buffer_r)
188  var=var_tmp(1)
189 
190  END SUBROUTINE bcast_omp_r
191 
192 
193  SUBROUTINE bcast_omp_r1(var)
194  IMPLICIT NONE
195  REAL,INTENT(INOUT) :: Var(:)
196 
197  CALL check_buffer_r(size(var))
198  CALL bcast_omp_rgen(var,size(var),buffer_r)
199 
200  END SUBROUTINE bcast_omp_r1
201 
202 
203  SUBROUTINE bcast_omp_r2(var)
204  IMPLICIT NONE
205  REAL,INTENT(INOUT) :: Var(:,:)
206 
207  CALL check_buffer_r(size(var))
208  CALL bcast_omp_rgen(var,size(var),buffer_r)
209 
210  END SUBROUTINE bcast_omp_r2
211 
212 
213  SUBROUTINE bcast_omp_r3(var)
214  IMPLICIT NONE
215  REAL,INTENT(INOUT) :: Var(:,:,:)
216 
217  CALL check_buffer_r(size(var))
218  CALL bcast_omp_rgen(var,size(var),buffer_r)
219 
220  END SUBROUTINE bcast_omp_r3
221 
222 
223  SUBROUTINE bcast_omp_r4(var)
224  IMPLICIT NONE
225  REAL,INTENT(INOUT) :: Var(:,:,:,:)
226 
227  CALL check_buffer_r(size(var))
228  CALL bcast_omp_rgen(var,size(var),buffer_r)
229 
230  END SUBROUTINE bcast_omp_r4
231 
232 
233 !! -- Les booleans -- !!
234 
235  SUBROUTINE bcast_omp_l(var)
236  IMPLICIT NONE
237  LOGICAL,INTENT(INOUT) :: Var
238  LOGICAL :: Var_tmp(1)
239 
240  var_tmp(1)=var
241  CALL check_buffer_l(1)
242  CALL bcast_omp_lgen(var_tmp,1,buffer_l)
243  var=var_tmp(1)
244 
245  END SUBROUTINE bcast_omp_l
246 
247 
248  SUBROUTINE bcast_omp_l1(var)
249  IMPLICIT NONE
250  LOGICAL,INTENT(INOUT) :: Var(:)
251 
252  CALL check_buffer_l(size(var))
253  CALL bcast_omp_lgen(var,size(var),buffer_l)
254 
255  END SUBROUTINE bcast_omp_l1
256 
257 
258  SUBROUTINE bcast_omp_l2(var)
259  IMPLICIT NONE
260  LOGICAL,INTENT(INOUT) :: Var(:,:)
261 
262  CALL check_buffer_l(size(var))
263  CALL bcast_omp_lgen(var,size(var),buffer_l)
264 
265  END SUBROUTINE bcast_omp_l2
266 
267 
268  SUBROUTINE bcast_omp_l3(var)
269  IMPLICIT NONE
270  LOGICAL,INTENT(INOUT) :: Var(:,:,:)
271 
272  CALL check_buffer_l(size(var))
273  CALL bcast_omp_lgen(var,size(var),buffer_l)
274 
275  END SUBROUTINE bcast_omp_l3
276 
277 
278  SUBROUTINE bcast_omp_l4(var)
279  IMPLICIT NONE
280  LOGICAL,INTENT(INOUT) :: Var(:,:,:,:)
281 
282  CALL check_buffer_l(size(var))
283  CALL bcast_omp_lgen(var,size(var),buffer_l)
284 
285  END SUBROUTINE bcast_omp_l4
286 
287 
288 
289 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
290 !! Definition des Scatter --> 4D !!
291 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
292 
293  SUBROUTINE scatter_omp_i(VarIn, VarOut)
294  IMPLICIT NONE
295 
296  INTEGER,INTENT(IN),DIMENSION(:) :: VarIn
297  INTEGER,INTENT(OUT),DIMENSION(:) :: VarOut
298 
299  CALL check_buffer_i(size(varin))
300  CALL scatter_omp_igen(varin,varout,1,buffer_i)
301 
302  END SUBROUTINE scatter_omp_i
303 
304 
305  SUBROUTINE scatter_omp_i1(VarIn, VarOut)
306  IMPLICIT NONE
307 
308  INTEGER,INTENT(IN),DIMENSION(:,:) :: VarIn
309  INTEGER,INTENT(OUT),DIMENSION(:,:) :: VarOut
310 
311  CALL check_buffer_i(size(varin))
312  CALL scatter_omp_igen(varin,varout,Size(varout,2),buffer_i)
313 
314  END SUBROUTINE scatter_omp_i1
315 
316 
317  SUBROUTINE scatter_omp_i2(VarIn, VarOut)
318  IMPLICIT NONE
319 
320  INTEGER,INTENT(IN),DIMENSION(:,:,:) :: VarIn
321  INTEGER,INTENT(OUT),DIMENSION(:,:,:) :: VarOut
322 
323  CALL check_buffer_i(size(varin))
324  CALL scatter_omp_igen(varin,varout,Size(varout,2)*Size(varout,3),buffer_i)
325 
326  END SUBROUTINE scatter_omp_i2
327 
328 
329  SUBROUTINE scatter_omp_i3(VarIn, VarOut)
330  IMPLICIT NONE
331 
332  INTEGER,INTENT(IN),DIMENSION(:,:,:,:) :: VarIn
333  INTEGER,INTENT(OUT),DIMENSION(:,:,:,:) :: VarOut
334 
335  CALL check_buffer_i(size(varin))
336  CALL scatter_omp_igen(varin,varout,Size(varout,2)*Size(varout,3)*Size(varout,4),buffer_i)
337 
338  END SUBROUTINE scatter_omp_i3
339 
340 
341 
342 
343  SUBROUTINE scatter_omp_r(VarIn, VarOut)
344  IMPLICIT NONE
345 
346  REAL,INTENT(IN),DIMENSION(:) :: VarIn
347  REAL,INTENT(OUT),DIMENSION(:) :: VarOut
348 
349  CALL check_buffer_r(size(varin))
350  CALL scatter_omp_rgen(varin,varout,1,buffer_r)
351 
352  END SUBROUTINE scatter_omp_r
353 
354 
355  SUBROUTINE scatter_omp_r1(VarIn, VarOut)
356  IMPLICIT NONE
357 
358  REAL,INTENT(IN),DIMENSION(:,:) :: VarIn
359  REAL,INTENT(OUT),DIMENSION(:,:) :: VarOut
360 
361  CALL check_buffer_r(size(varin))
362  CALL scatter_omp_rgen(varin,varout,Size(varout,2),buffer_r)
363 
364  END SUBROUTINE scatter_omp_r1
365 
366 
367  SUBROUTINE scatter_omp_r2(VarIn, VarOut)
368  IMPLICIT NONE
369 
370  REAL,INTENT(IN),DIMENSION(:,:,:) :: VarIn
371  REAL,INTENT(OUT),DIMENSION(:,:,:) :: VarOut
372 
373  CALL check_buffer_r(size(varin))
374  CALL scatter_omp_rgen(varin,varout,Size(varout,2)*Size(varout,3),buffer_r)
375 
376  END SUBROUTINE scatter_omp_r2
377 
378 
379  SUBROUTINE scatter_omp_r3(VarIn, VarOut)
380  IMPLICIT NONE
381 
382  REAL,INTENT(IN),DIMENSION(:,:,:,:) :: VarIn
383  REAL,INTENT(OUT),DIMENSION(:,:,:,:) :: VarOut
384 
385  CALL check_buffer_r(size(varin))
386  CALL scatter_omp_rgen(varin,varout,Size(varout,2)*Size(varout,3)*Size(varout,4),buffer_r)
387 
388  END SUBROUTINE scatter_omp_r3
389 
390 
391 
392  SUBROUTINE scatter_omp_l(VarIn, VarOut)
393  IMPLICIT NONE
394 
395  LOGICAL,INTENT(IN),DIMENSION(:) :: VarIn
396  LOGICAL,INTENT(OUT),DIMENSION(:) :: VarOut
397 
398  CALL check_buffer_l(size(varin))
399  CALL scatter_omp_lgen(varin,varout,1,buffer_l)
400 
401  END SUBROUTINE scatter_omp_l
402 
403 
404  SUBROUTINE scatter_omp_l1(VarIn, VarOut)
405  IMPLICIT NONE
406 
407  LOGICAL,INTENT(IN),DIMENSION(:,:) :: VarIn
408  LOGICAL,INTENT(OUT),DIMENSION(:,:) :: VarOut
409 
410  CALL check_buffer_l(size(varin))
411  CALL scatter_omp_lgen(varin,varout,Size(varout,2),buffer_l)
412 
413  END SUBROUTINE scatter_omp_l1
414 
415 
416  SUBROUTINE scatter_omp_l2(VarIn, VarOut)
417  IMPLICIT NONE
418 
419  LOGICAL,INTENT(IN),DIMENSION(:,:,:) :: VarIn
420  LOGICAL,INTENT(OUT),DIMENSION(:,:,:) :: VarOut
421 
422  CALL check_buffer_l(size(varin))
423  CALL scatter_omp_lgen(varin,varout,Size(varout,2)*Size(varout,3),buffer_l)
424 
425  END SUBROUTINE scatter_omp_l2
426 
427 
428  SUBROUTINE scatter_omp_l3(VarIn, VarOut)
429  IMPLICIT NONE
430 
431  LOGICAL,INTENT(IN),DIMENSION(:,:,:,:) :: VarIn
432  LOGICAL,INTENT(OUT),DIMENSION(:,:,:,:) :: VarOut
433 
434  CALL check_buffer_l(size(varin))
435  CALL scatter_omp_lgen(varin,varout,Size(varout,2)*Size(varout,3)*Size(varout,4),buffer_l)
436 
437  END SUBROUTINE scatter_omp_l3
438 
439 
440  SUBROUTINE gather_omp_i(VarIn, VarOut)
441  IMPLICIT NONE
442 
443  INTEGER,INTENT(IN),DIMENSION(:) :: VarIn
444  INTEGER,INTENT(OUT),DIMENSION(:) :: VarOut
445 
446  CALL check_buffer_i(size(varout))
447  CALL gather_omp_igen(varin,varout,1,buffer_i)
448 
449  END SUBROUTINE gather_omp_i
450 
451 
452  SUBROUTINE gather_omp_i1(VarIn, VarOut)
453  IMPLICIT NONE
454 
455  INTEGER,INTENT(IN),DIMENSION(:,:) :: VarIn
456  INTEGER,INTENT(OUT),DIMENSION(:,:) :: VarOut
457 
458  CALL check_buffer_i(size(varout))
459  CALL gather_omp_igen(varin,varout,Size(varin,2),buffer_i)
460 
461  END SUBROUTINE gather_omp_i1
462 
463 
464  SUBROUTINE gather_omp_i2(VarIn, VarOut)
465  IMPLICIT NONE
466 
467  INTEGER,INTENT(IN),DIMENSION(:,:,:) :: VarIn
468  INTEGER,INTENT(OUT),DIMENSION(:,:,:) :: VarOut
469 
470  CALL check_buffer_i(size(varout))
471  CALL gather_omp_igen(varin,varout,Size(varin,2)*Size(varin,3),buffer_i)
472 
473  END SUBROUTINE gather_omp_i2
474 
475 
476  SUBROUTINE gather_omp_i3(VarIn, VarOut)
477  IMPLICIT NONE
478 
479  INTEGER,INTENT(IN),DIMENSION(:,:,:,:) :: VarIn
480  INTEGER,INTENT(OUT),DIMENSION(:,:,:,:) :: VarOut
481 
482  CALL check_buffer_i(size(varout))
483  CALL gather_omp_igen(varin,varout,Size(varin,2)*Size(varin,3)*Size(varin,4),buffer_i)
484 
485  END SUBROUTINE gather_omp_i3
486 
487 
488 
489  SUBROUTINE gather_omp_r(VarIn, VarOut)
490  IMPLICIT NONE
491 
492  REAL,INTENT(IN),DIMENSION(:) :: VarIn
493  REAL,INTENT(OUT),DIMENSION(:) :: VarOut
494 
495  CALL check_buffer_r(size(varout))
496  CALL gather_omp_rgen(varin,varout,1,buffer_r)
497 
498  END SUBROUTINE gather_omp_r
499 
500 
501  SUBROUTINE gather_omp_r1(VarIn, VarOut)
502  IMPLICIT NONE
503 
504  REAL,INTENT(IN),DIMENSION(:,:) :: VarIn
505  REAL,INTENT(OUT),DIMENSION(:,:) :: VarOut
506 
507  CALL check_buffer_r(size(varout))
508  CALL gather_omp_rgen(varin,varout,Size(varin,2),buffer_r)
509 
510  END SUBROUTINE gather_omp_r1
511 
512 
513  SUBROUTINE gather_omp_r2(VarIn, VarOut)
514  IMPLICIT NONE
515 
516  REAL,INTENT(IN),DIMENSION(:,:,:) :: VarIn
517  REAL,INTENT(OUT),DIMENSION(:,:,:) :: VarOut
518 
519  CALL check_buffer_r(size(varout))
520  CALL gather_omp_rgen(varin,varout,Size(varin,2)*Size(varin,3),buffer_r)
521 
522  END SUBROUTINE gather_omp_r2
523 
524 
525  SUBROUTINE gather_omp_r3(VarIn, VarOut)
526  IMPLICIT NONE
527 
528  REAL,INTENT(IN),DIMENSION(:,:,:,:) :: VarIn
529  REAL,INTENT(OUT),DIMENSION(:,:,:,:) :: VarOut
530 
531  CALL check_buffer_r(size(varout))
532  CALL gather_omp_rgen(varin,varout,Size(varin,2)*Size(varin,3)*Size(varin,4),buffer_r)
533 
534  END SUBROUTINE gather_omp_r3
535 
536 
537  SUBROUTINE gather_omp_l(VarIn, VarOut)
538  IMPLICIT NONE
539 
540  LOGICAL,INTENT(IN),DIMENSION(:) :: VarIn
541  LOGICAL,INTENT(OUT),DIMENSION(:) :: VarOut
542 
543  CALL check_buffer_l(size(varout))
544  CALL gather_omp_lgen(varin,varout,1,buffer_l)
545 
546  END SUBROUTINE gather_omp_l
547 
548 
549  SUBROUTINE gather_omp_l1(VarIn, VarOut)
550  IMPLICIT NONE
551 
552  LOGICAL,INTENT(IN),DIMENSION(:,:) :: VarIn
553  LOGICAL,INTENT(OUT),DIMENSION(:,:) :: VarOut
554 
555  CALL check_buffer_l(size(varout))
556  CALL gather_omp_lgen(varin,varout,Size(varin,2),buffer_l)
557 
558  END SUBROUTINE gather_omp_l1
559 
560 
561  SUBROUTINE gather_omp_l2(VarIn, VarOut)
562  IMPLICIT NONE
563 
564  LOGICAL,INTENT(IN),DIMENSION(:,:,:) :: VarIn
565  LOGICAL,INTENT(OUT),DIMENSION(:,:,:) :: VarOut
566 
567  CALL check_buffer_l(size(varout))
568  CALL gather_omp_lgen(varin,varout,Size(varin,2)*Size(varin,3),buffer_l)
569 
570  END SUBROUTINE gather_omp_l2
571 
572 
573  SUBROUTINE gather_omp_l3(VarIn, VarOut)
574  IMPLICIT NONE
575 
576  LOGICAL,INTENT(IN),DIMENSION(:,:,:,:) :: VarIn
577  LOGICAL,INTENT(OUT),DIMENSION(:,:,:,:) :: VarOut
578 
579  CALL check_buffer_l(size(varout))
580  CALL gather_omp_lgen(varin,varout,Size(varin,2)*Size(varin,3)*Size(varin,4),buffer_l)
581 
582  END SUBROUTINE gather_omp_l3
583 
584 
585 
586 
587  SUBROUTINE reduce_sum_omp_i(VarIn, VarOut)
588  IMPLICIT NONE
589 
590  INTEGER,INTENT(IN) :: VarIn
591  INTEGER,INTENT(OUT) :: VarOut
592  INTEGER :: VarIn_tmp(1)
593  INTEGER :: VarOut_tmp(1)
594 
595  varin_tmp(1)=varin
596  CALL check_buffer_i(1)
597  CALL reduce_sum_omp_igen(varin_tmp,varout_tmp,1,buffer_i)
598  varout=varout_tmp(1)
599 
600  END SUBROUTINE reduce_sum_omp_i
601 
602  SUBROUTINE reduce_sum_omp_i1(VarIn, VarOut)
603  IMPLICIT NONE
604 
605  INTEGER,INTENT(IN),DIMENSION(:) :: VarIn
606  INTEGER,INTENT(OUT),DIMENSION(:) :: VarOut
607 
608  CALL check_buffer_i(size(varin))
609  CALL reduce_sum_omp_igen(varin,varout,Size(varin),buffer_i)
610 
611  END SUBROUTINE reduce_sum_omp_i1
612 
613 
614  SUBROUTINE reduce_sum_omp_i2(VarIn, VarOut)
615  IMPLICIT NONE
616 
617  INTEGER,INTENT(IN),DIMENSION(:,:) :: VarIn
618  INTEGER,INTENT(OUT),DIMENSION(:,:) :: VarOut
619 
620  CALL check_buffer_i(size(varin))
621  CALL reduce_sum_omp_igen(varin,varout,Size(varin),buffer_i)
622 
623  END SUBROUTINE reduce_sum_omp_i2
624 
625 
626  SUBROUTINE reduce_sum_omp_i3(VarIn, VarOut)
627  IMPLICIT NONE
628 
629  INTEGER,INTENT(IN),DIMENSION(:,:,:) :: VarIn
630  INTEGER,INTENT(OUT),DIMENSION(:,:,:) :: VarOut
631 
632  CALL check_buffer_i(size(varin))
633  CALL reduce_sum_omp_igen(varin,varout,Size(varin),buffer_i)
634 
635  END SUBROUTINE reduce_sum_omp_i3
636 
637 
638  SUBROUTINE reduce_sum_omp_i4(VarIn, VarOut)
639  IMPLICIT NONE
640 
641  INTEGER,INTENT(IN),DIMENSION(:,:,:,:) :: VarIn
642  INTEGER,INTENT(OUT),DIMENSION(:,:,:,:) :: VarOut
643 
644  CALL check_buffer_i(size(varin))
645  CALL reduce_sum_omp_igen(varin,varout,Size(varin),buffer_i)
646 
647  END SUBROUTINE reduce_sum_omp_i4
648 
649 
650  SUBROUTINE reduce_sum_omp_r(VarIn, VarOut)
651  IMPLICIT NONE
652 
653  REAL,INTENT(IN) :: VarIn
654  REAL,INTENT(OUT) :: VarOut
655  REAL :: VarIn_tmp(1)
656  REAL :: VarOut_tmp(1)
657 
658  varin_tmp(1)=varin
659  CALL check_buffer_r(1)
660  CALL reduce_sum_omp_rgen(varin_tmp,varout_tmp,1,buffer_r)
661  varout=varout_tmp(1)
662 
663  END SUBROUTINE reduce_sum_omp_r
664 
665  SUBROUTINE reduce_sum_omp_r1(VarIn, VarOut)
666  IMPLICIT NONE
667 
668  REAL,INTENT(IN),DIMENSION(:) :: VarIn
669  REAL,INTENT(OUT),DIMENSION(:) :: VarOut
670 
671  CALL check_buffer_r(size(varin))
672  CALL reduce_sum_omp_rgen(varin,varout,Size(varin),buffer_r)
673 
674  END SUBROUTINE reduce_sum_omp_r1
675 
676 
677  SUBROUTINE reduce_sum_omp_r2(VarIn, VarOut)
678  IMPLICIT NONE
679 
680  REAL,INTENT(IN),DIMENSION(:,:) :: VarIn
681  REAL,INTENT(OUT),DIMENSION(:,:) :: VarOut
682 
683  CALL check_buffer_r(size(varin))
684  CALL reduce_sum_omp_rgen(varin,varout,Size(varin),buffer_r)
685 
686  END SUBROUTINE reduce_sum_omp_r2
687 
688 
689  SUBROUTINE reduce_sum_omp_r3(VarIn, VarOut)
690  IMPLICIT NONE
691 
692  REAL,INTENT(IN),DIMENSION(:,:,:) :: VarIn
693  REAL,INTENT(OUT),DIMENSION(:,:,:) :: VarOut
694 
695  CALL check_buffer_r(size(varin))
696  CALL reduce_sum_omp_rgen(varin,varout,Size(varin),buffer_r)
697 
698  END SUBROUTINE reduce_sum_omp_r3
699 
700 
701  SUBROUTINE reduce_sum_omp_r4(VarIn, VarOut)
702  IMPLICIT NONE
703 
704  REAL,INTENT(IN),DIMENSION(:,:,:,:) :: VarIn
705  REAL,INTENT(OUT),DIMENSION(:,:,:,:) :: VarOut
706 
707  CALL check_buffer_r(size(varin))
708  CALL reduce_sum_omp_rgen(varin,varout,Size(varin),buffer_r)
709 
710  END SUBROUTINE reduce_sum_omp_r4
711 
712 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
713 ! LES ROUTINES GENERIQUES !
714 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
715 
716  SUBROUTINE bcast_omp_cgen(Var,Nb,Buff)
717  IMPLICIT NONE
718 
719  CHARACTER(LEN=*),INTENT(INOUT) :: Var
720  CHARACTER(LEN=*),INTENT(INOUT) :: Buff
721  INTEGER,INTENT(IN) :: Nb
722 
723  INTEGER :: i
724 
725  !$OMP MASTER
726  buff=var
727  !$OMP END MASTER
728  !$OMP BARRIER
729 
730  DO i=1,nb
731  var=buff
732  ENDDO
733  !$OMP BARRIER
734 
735  END SUBROUTINE bcast_omp_cgen
736 
737 
738 
739  SUBROUTINE bcast_omp_igen(Var,Nb,Buff)
740  IMPLICIT NONE
741 
742  INTEGER,INTENT(IN) :: Nb
743  INTEGER,DIMENSION(Nb),INTENT(INOUT) :: Var
744  INTEGER,DIMENSION(Nb),INTENT(INOUT) :: Buff
745 
746  INTEGER :: i
747 
748  !$OMP MASTER
749  DO i=1,nb
750  buff(i)=var(i)
751  ENDDO
752  !$OMP END MASTER
753  !$OMP BARRIER
754 
755  DO i=1,nb
756  var(i)=buff(i)
757  ENDDO
758  !$OMP BARRIER
759 
760  END SUBROUTINE bcast_omp_igen
761 
762 
763  SUBROUTINE bcast_omp_rgen(Var,Nb,Buff)
764  IMPLICIT NONE
765 
766  INTEGER,INTENT(IN) :: Nb
767  REAL,DIMENSION(Nb),INTENT(INOUT) :: Var
768  REAL,DIMENSION(Nb),INTENT(INOUT) :: Buff
769 
770  INTEGER :: i
771 
772  !$OMP MASTER
773  DO i=1,nb
774  buff(i)=var(i)
775  ENDDO
776  !$OMP END MASTER
777  !$OMP BARRIER
778 
779  DO i=1,nb
780  var(i)=buff(i)
781  ENDDO
782  !$OMP BARRIER
783 
784  END SUBROUTINE bcast_omp_rgen
785 
786  SUBROUTINE bcast_omp_lgen(Var,Nb,Buff)
787  IMPLICIT NONE
788 
789  INTEGER,INTENT(IN) :: Nb
790  LOGICAL,DIMENSION(Nb),INTENT(INOUT) :: Var
791  LOGICAL,DIMENSION(Nb),INTENT(INOUT) :: Buff
792 
793  INTEGER :: i
794 
795  !$OMP MASTER
796  DO i=1,nb
797  buff(i)=var(i)
798  ENDDO
799  !$OMP END MASTER
800  !$OMP BARRIER
801 
802  DO i=1,nb
803  var(i)=buff(i)
804  ENDDO
805  !$OMP BARRIER
806 
807  END SUBROUTINE bcast_omp_lgen
808 
809 
810  SUBROUTINE scatter_omp_igen(VarIn,VarOut,dimsize,Buff)
812  USE mod_phys_lmdz_mpi_data, ONLY : klon_mpi
813  IMPLICIT NONE
814 
815  INTEGER,INTENT(IN) :: dimsize
816  INTEGER,INTENT(IN),DIMENSION(klon_mpi,dimsize) :: VarIn
817  INTEGER,INTENT(OUT),DIMENSION(klon_omp,dimsize) :: VarOut
818  INTEGER,INTENT(INOUT),DIMENSION(klon_mpi,dimsize) :: Buff
819 
820  INTEGER :: i,ij
821 
822  !$OMP MASTER
823  DO i=1,dimsize
824  DO ij=1,klon_mpi
825  buff(ij,i)=varin(ij,i)
826  ENDDO
827  ENDDO
828  !$OMP END MASTER
829  !$OMP BARRIER
830 
831  DO i=1,dimsize
832  DO ij=1,klon_omp
833  varout(ij,i)=buff(klon_omp_begin-1+ij,i)
834  ENDDO
835  ENDDO
836  !$OMP BARRIER
837 
838  END SUBROUTINE scatter_omp_igen
839 
840 
841  SUBROUTINE scatter_omp_rgen(VarIn,VarOut,dimsize,Buff)
843  USE mod_phys_lmdz_mpi_data, ONLY : klon_mpi
844  IMPLICIT NONE
845 
846  INTEGER,INTENT(IN) :: dimsize
847  REAL,INTENT(IN),DIMENSION(klon_mpi,dimsize) :: VarIn
848  REAL,INTENT(OUT),DIMENSION(klon_omp,dimsize) :: VarOut
849  REAL,INTENT(INOUT),DIMENSION(klon_mpi,dimsize) :: Buff
850 
851  INTEGER :: i,ij
852 
853  !$OMP MASTER
854  DO i=1,dimsize
855  DO ij=1,klon_mpi
856  buff(ij,i)=varin(ij,i)
857  ENDDO
858  ENDDO
859  !$OMP END MASTER
860  !$OMP BARRIER
861 
862  DO i=1,dimsize
863  DO ij=1,klon_omp
864  varout(ij,i)=buff(klon_omp_begin-1+ij,i)
865  ENDDO
866  ENDDO
867  !$OMP BARRIER
868 
869  END SUBROUTINE scatter_omp_rgen
870 
871 
872  SUBROUTINE scatter_omp_lgen(VarIn,VarOut,dimsize,Buff)
874  USE mod_phys_lmdz_mpi_data, ONLY : klon_mpi
875  IMPLICIT NONE
876 
877  INTEGER,INTENT(IN) :: dimsize
878  LOGICAL,INTENT(IN),DIMENSION(klon_mpi,dimsize) :: VarIn
879  LOGICAL,INTENT(OUT),DIMENSION(klon_omp,dimsize) :: VarOut
880  LOGICAL,INTENT(INOUT),DIMENSION(klon_mpi,dimsize) :: Buff
881 
882  INTEGER :: i,ij
883 
884  !$OMP MASTER
885  DO i=1,dimsize
886  DO ij=1,klon_mpi
887  buff(ij,i)=varin(ij,i)
888  ENDDO
889  ENDDO
890  !$OMP END MASTER
891  !$OMP BARRIER
892 
893  DO i=1,dimsize
894  DO ij=1,klon_omp
895  varout(ij,i)=buff(klon_omp_begin-1+ij,i)
896  ENDDO
897  ENDDO
898  !$OMP BARRIER
899 
900  END SUBROUTINE scatter_omp_lgen
901 
902 
903 
904 
905 
906  SUBROUTINE gather_omp_igen(VarIn,VarOut,dimsize,Buff)
908  USE mod_phys_lmdz_mpi_data, ONLY : klon_mpi
909  IMPLICIT NONE
910 
911  INTEGER,INTENT(IN) :: dimsize
912  INTEGER,INTENT(IN),DIMENSION(klon_omp,dimsize) :: VarIn
913  INTEGER,INTENT(OUT),DIMENSION(klon_mpi,dimsize) :: VarOut
914  INTEGER,INTENT(INOUT),DIMENSION(klon_mpi,dimsize) :: Buff
915 
916  INTEGER :: i,ij
917 
918  DO i=1,dimsize
919  DO ij=1,klon_omp
920  buff(klon_omp_begin-1+ij,i)=varin(ij,i)
921  ENDDO
922  ENDDO
923  !$OMP BARRIER
924 
925 
926  !$OMP MASTER
927  DO i=1,dimsize
928  DO ij=1,klon_mpi
929  varout(ij,i)=buff(ij,i)
930  ENDDO
931  ENDDO
932  !$OMP END MASTER
933  !$OMP BARRIER
934 
935  END SUBROUTINE gather_omp_igen
936 
937 
938  SUBROUTINE gather_omp_rgen(VarIn,VarOut,dimsize,Buff)
940  USE mod_phys_lmdz_mpi_data, ONLY : klon_mpi
941  IMPLICIT NONE
942 
943  INTEGER,INTENT(IN) :: dimsize
944  REAL,INTENT(IN),DIMENSION(klon_omp,dimsize) :: VarIn
945  REAL,INTENT(OUT),DIMENSION(klon_mpi,dimsize) :: VarOut
946  REAL,INTENT(INOUT),DIMENSION(klon_mpi,dimsize) :: Buff
947 
948  INTEGER :: i,ij
949 
950  DO i=1,dimsize
951  DO ij=1,klon_omp
952  buff(klon_omp_begin-1+ij,i)=varin(ij,i)
953  ENDDO
954  ENDDO
955  !$OMP BARRIER
956 
957 
958  !$OMP MASTER
959  DO i=1,dimsize
960  DO ij=1,klon_mpi
961  varout(ij,i)=buff(ij,i)
962  ENDDO
963  ENDDO
964  !$OMP END MASTER
965  !$OMP BARRIER
966 
967  END SUBROUTINE gather_omp_rgen
968 
969 
970  SUBROUTINE gather_omp_lgen(VarIn,VarOut,dimsize,Buff)
972  USE mod_phys_lmdz_mpi_data, ONLY : klon_mpi
973  IMPLICIT NONE
974 
975  INTEGER,INTENT(IN) :: dimsize
976  LOGICAL,INTENT(IN),DIMENSION(klon_omp,dimsize) :: VarIn
977  LOGICAL,INTENT(OUT),DIMENSION(klon_mpi,dimsize) :: VarOut
978  LOGICAL,INTENT(INOUT),DIMENSION(klon_mpi,dimsize) :: Buff
979 
980  INTEGER :: i,ij
981 
982  DO i=1,dimsize
983  DO ij=1,klon_omp
984  buff(klon_omp_begin-1+ij,i)=varin(ij,i)
985  ENDDO
986  ENDDO
987  !$OMP BARRIER
988 
989 
990  !$OMP MASTER
991  DO i=1,dimsize
992  DO ij=1,klon_mpi
993  varout(ij,i)=buff(ij,i)
994  ENDDO
995  ENDDO
996  !$OMP END MASTER
997  !$OMP BARRIER
998 
999  END SUBROUTINE gather_omp_lgen
1000 
1001 
1002  SUBROUTINE reduce_sum_omp_igen(VarIn,VarOut,dimsize,Buff)
1003  IMPLICIT NONE
1004 
1005  INTEGER,INTENT(IN) :: dimsize
1006  INTEGER,INTENT(IN),DIMENSION(dimsize) :: VarIn
1007  INTEGER,INTENT(OUT),DIMENSION(dimsize) :: VarOut
1008  INTEGER,INTENT(INOUT),DIMENSION(dimsize) :: Buff
1009 
1010  INTEGER :: i
1011 
1012  !$OMP MASTER
1013  buff(:)=0
1014  !$OMP END MASTER
1015  !$OMP BARRIER
1016 
1017  !$OMP CRITICAL
1018  DO i=1,dimsize
1019  buff(i)=buff(i)+varin(i)
1020  ENDDO
1021  !$OMP END CRITICAL
1022  !$OMP BARRIER
1023 
1024  !$OMP MASTER
1025  DO i=1,dimsize
1026  varout(i)=buff(i)
1027  ENDDO
1028  !$OMP END MASTER
1029  !$OMP BARRIER
1030 
1031  END SUBROUTINE reduce_sum_omp_igen
1032 
1033  SUBROUTINE reduce_sum_omp_rgen(VarIn,VarOut,dimsize,Buff)
1034  IMPLICIT NONE
1035 
1036  INTEGER,INTENT(IN) :: dimsize
1037  REAL,INTENT(IN),DIMENSION(dimsize) :: VarIn
1038  REAL,INTENT(OUT),DIMENSION(dimsize) :: VarOut
1039  REAL,INTENT(INOUT),DIMENSION(dimsize) :: Buff
1040 
1041  INTEGER :: i
1042 
1043  !$OMP MASTER
1044  buff(:)=0
1045  !$OMP END MASTER
1046  !$OMP BARRIER
1047 
1048  !$OMP CRITICAL
1049  DO i=1,dimsize
1050  buff(i)=buff(i)+varin(i)
1051  ENDDO
1052  !$OMP END CRITICAL
1053  !$OMP BARRIER
1054 
1055  !$OMP MASTER
1056  DO i=1,dimsize
1057  varout(i)=buff(i)
1058  ENDDO
1059  !$OMP END MASTER
1060  !$OMP BARRIER
1061 
1062  END SUBROUTINE reduce_sum_omp_rgen
1063 
1064 END MODULE mod_phys_lmdz_omp_transfert
integer, dimension(:), allocatable, save buffer_i
subroutine scatter_omp_l1(VarIn, VarOut)
subroutine reduce_sum_omp_r4(VarIn, VarOut)
subroutine bcast_omp_lgen(Var, Nb, Buff)
subroutine bcast_omp_cgen(Var, Nb, Buff)
subroutine bcast_omp_igen(Var, Nb, Buff)
subroutine scatter_omp_i3(VarIn, VarOut)
subroutine scatter_omp_lgen(VarIn, VarOut, dimsize, Buff)
logical, dimension(:), allocatable, save buffer_l
subroutine reduce_sum_omp_r2(VarIn, VarOut)
subroutine reduce_sum_omp_i1(VarIn, VarOut)
subroutine reduce_sum_omp_igen(VarIn, VarOut, dimsize, Buff)
subroutine scatter_omp_r3(VarIn, VarOut)
subroutine gather_omp_rgen(VarIn, VarOut, dimsize, Buff)
subroutine reduce_sum_omp_i3(VarIn, VarOut)
subroutine scatter_omp_i2(VarIn, VarOut)
subroutine scatter_omp_r1(VarIn, VarOut)
subroutine reduce_sum_omp_i(VarIn, VarOut)
subroutine reduce_sum_omp_r(VarIn, VarOut)
subroutine scatter_omp_r2(VarIn, VarOut)
subroutine scatter_omp_l3(VarIn, VarOut)
subroutine gather_omp_igen(VarIn, VarOut, dimsize, Buff)
subroutine reduce_sum_omp_i4(VarIn, VarOut)
subroutine reduce_sum_omp_r1(VarIn, VarOut)
subroutine reduce_sum_omp_rgen(VarIn, VarOut, dimsize, Buff)
subroutine scatter_omp_igen(VarIn, VarOut, dimsize, Buff)
subroutine scatter_omp_l2(VarIn, VarOut)
subroutine gather_omp_lgen(VarIn, VarOut, dimsize, Buff)
subroutine bcast_omp_rgen(Var, Nb, Buff)
subroutine scatter_omp_i1(VarIn, VarOut)
subroutine reduce_sum_omp_i2(VarIn, VarOut)
subroutine scatter_omp_rgen(VarIn, VarOut, dimsize, Buff)
real, dimension(:), allocatable, save buffer_r
character(len=size_min), save buffer_c
subroutine reduce_sum_omp_r3(VarIn, VarOut)