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