source: codes/icosagcm/trunk/src/dissip_gcm.f90 @ 151

Last change on this file since 151 was 151, checked in by ymipsl, 11 years ago

Implementation of mixte parallelism MPI/OpenMP into src directory

YM

File size: 24.2 KB
Line 
1MODULE dissip_gcm_mod
2  USE icosa
3
4  PRIVATE
5
6  TYPE(t_field),POINTER,SAVE :: f_due_diss1(:)
7  TYPE(t_field),POINTER,SAVE :: f_due_diss2(:)
8
9  TYPE(t_field),POINTER,SAVE :: f_theta(:), f_phi(:), f_pk(:), f_pks(:), f_p(:)
10  TYPE(t_field),POINTER,SAVE :: f_dtheta_diss(:)
11  TYPE(t_field),POINTER,SAVE :: f_dtheta_rhodz_diss(:)
12  TYPE(t_message) :: req_due, req_dtheta 
13 
14  INTEGER,SAVE :: nitergdiv=1
15  INTEGER,SAVE :: nitergrot=1
16  INTEGER,SAVE :: niterdivgrad=1
17
18  REAL,ALLOCATABLE,SAVE :: tau_graddiv(:)
19  REAL,ALLOCATABLE,SAVE :: tau_gradrot(:)
20  REAL,ALLOCATABLE,SAVE :: tau_divgrad(:)
21 
22  REAL,SAVE :: cgraddiv
23  REAL,SAVE :: cgradrot
24  REAL,SAVE :: cdivgrad
25
26  INTEGER, SAVE :: rayleigh_friction_type, rayleigh_shear
27  REAL, save    :: rayleigh_tau
28 
29!  INTEGER,SAVE :: itau_dissip
30  REAL,SAVE    :: dtdissip
31 
32  PUBLIC init_dissip, dissip
33CONTAINS
34
35  SUBROUTINE allocate_dissip
36  USE icosa
37  IMPLICIT NONE 
38    CALL allocate_field(f_due_diss1,field_u,type_real,llm)
39    CALL allocate_field(f_due_diss2,field_u,type_real,llm)
40    CALL allocate_field(f_theta,field_t,type_real,llm)
41    CALL allocate_field(f_dtheta_diss,field_t,type_real,llm)
42    CALL allocate_field(f_dtheta_rhodz_diss,field_t,type_real,llm)
43
44    CALL allocate_field(f_phi,field_t,type_real,llm)
45    CALL allocate_field(f_pk,field_t,type_real,llm)
46    CALL allocate_field(f_p,field_t,type_real,llm+1)
47    CALL allocate_field(f_pks,field_t,type_real)
48   
49    ALLOCATE(tau_graddiv(llm))
50    ALLOCATE(tau_gradrot(llm))   
51    ALLOCATE(tau_divgrad(llm))
52  END SUBROUTINE allocate_dissip
53 
54  SUBROUTINE init_dissip
55  USE icosa
56  USE disvert_mod
57  USE mpi_mod
58  USE mpipara
59  USE transfert_mod
60  USE time_mod
61 
62  IMPLICIT NONE
63 
64  TYPE(t_field),POINTER :: f_u(:)
65  TYPE(t_field),POINTER :: f_du(:)
66  REAL(rstd),POINTER    :: u(:)
67  REAL(rstd),POINTER    :: du(:)
68  TYPE(t_field),POINTER :: f_theta(:)
69  TYPE(t_field),POINTER :: f_dtheta(:)
70  REAL(rstd),POINTER    :: theta(:)
71  REAL(rstd),POINTER    :: dtheta(:)
72  REAL(rstd)            :: dumax,dumax1
73  REAL(rstd)            :: dthetamax,dthetamax1
74  REAL(rstd)            :: r
75  REAL(rstd)            :: tau
76  REAL(rstd)            :: zz, zvert, fact
77  INTEGER               :: l
78  CHARACTER(len=255)    :: rayleigh_friction_key
79  REAL(rstd)            :: mintau
80 
81           
82  INTEGER :: i,j,n,ind,it,iter
83
84   rayleigh_friction_key='none'
85   CALL getin("rayleigh_friction_type",rayleigh_friction_key)
86   SELECT CASE(TRIM(rayleigh_friction_key))
87   CASE('none')
88      rayleigh_friction_type=0
89      IF (is_mpi_root) PRINT *, 'No Rayleigh friction'
90   CASE('dcmip2_schaer_noshear')
91      rayleigh_friction_type=1
92      rayleigh_shear=0
93      IF (is_mpi_root) PRINT *, 'Rayleigh friction : Schaer-like mountain without shear DCMIP2.1'
94   CASE('dcmip2_schaer_shear')
95      rayleigh_shear=1
96      rayleigh_friction_type=2
97      IF (is_mpi_root) PRINT *, 'Rayleigh friction : Schaer-like mountain with shear DCMIP2.2'
98   CASE DEFAULT
99      IF (is_mpi_root) PRINT *, 'Bad selector : rayleigh_friction_type =', TRIM(rayleigh_friction_key), ' in dissip_gcm.f90/init_dissip'
100      STOP
101   END SELECT
102
103   IF(rayleigh_friction_type>0) THEN
104      rayleigh_tau=0.
105      CALL getin("rayleigh_friction_tau",rayleigh_tau)
106      rayleigh_tau = rayleigh_tau / scale_factor
107      IF(rayleigh_tau<=0) THEN
108         IF (is_mpi_root) PRINT *, 'Forbidden : negative value for rayleigh_friction_tau =',rayleigh_tau
109         STOP
110      END IF
111   END IF
112
113    CALL allocate_dissip
114    CALL allocate_field(f_u,field_u,type_real)
115    CALL allocate_field(f_du,field_u,type_real)
116    CALL allocate_field(f_theta,field_t,type_real)
117    CALL allocate_field(f_dtheta,field_t,type_real)
118   
119    CALL init_message(f_due_diss1,req_e1_vect,req_due)
120    CALL init_message(f_dtheta_diss,req_i1,req_dtheta)
121
122    tau_graddiv(:)=5000
123    CALL getin("tau_graddiv",tau)
124    tau_graddiv(:)=tau/scale_factor
125
126    CALL getin("nitergdiv",nitergdiv)
127   
128    tau_gradrot(:)=5000
129    CALL getin("tau_gradrot",tau_gradrot)
130    tau_gradrot(:)=tau/scale_factor
131
132    CALL getin("nitergrot",nitergrot)
133   
134
135    tau_divgrad(:)=5000
136    CALL getin("tau_divgrad",tau)
137    tau_divgrad(:)=tau/scale_factor
138
139    CALL getin("niterdivgrad",niterdivgrad)
140
141
142    cgraddiv=-1
143    cdivgrad=-1
144    cgradrot=-1
145 
146    CALL RANDOM_SEED()
147   
148    DO ind=1,ndomain
149      CALL swap_dimensions(ind)
150      CALL swap_geometry(ind)
151      u=f_u(ind)
152
153      DO j=jj_begin,jj_end
154        DO i=ii_begin,ii_end
155          n=(j-1)*iim+i
156          CALL RANDOM_NUMBER(r)
157          u(n+u_right)=r-0.5
158          CALL RANDOM_NUMBER(r)
159          u(n+u_lup)=r-0.5
160          CALL RANDOM_NUMBER(r)
161          u(n+u_ldown)=r-0.5
162       ENDDO
163      ENDDO
164       
165    ENDDO
166 
167
168
169    DO it=1,20
170     
171      dumax=0
172      DO iter=1,nitergdiv
173        CALL transfert_request(f_u,req_e1_vect)
174        DO ind=1,ndomain
175          CALL swap_dimensions(ind)
176          CALL swap_geometry(ind)
177          u=f_u(ind)
178          du=f_du(ind)
179          CALL compute_gradiv(u,du,1,1)
180          u=du
181        ENDDO
182      ENDDO
183
184      CALL transfert_request(f_du,req_e1_vect)
185     
186      DO ind=1,ndomain
187        CALL swap_dimensions(ind)
188        CALL swap_geometry(ind)
189        u=f_u(ind)
190        du=f_du(ind)
191         
192        DO j=jj_begin,jj_end
193          DO i=ii_begin,ii_end
194            n=(j-1)*iim+i
195            if (le(n+u_right)>1e-100) dumax=MAX(dumax,ABS(du(n+u_right)))
196            if (le(n+u_lup)>1e-100) dumax=MAX(dumax,ABS(du(n+u_lup)))
197            if (le(n+u_ldown)>1e-100) dumax=MAX(dumax,ABS(du(n+u_ldown)))
198          ENDDO
199        ENDDO
200      ENDDO
201
202      IF (using_mpi) THEN
203        CALL MPI_ALLREDUCE(dumax,dumax1,1,MPI_REAL8,MPI_MAX,comm_icosa,ierr)
204        dumax=dumax1
205      ENDIF
206                       
207      DO ind=1,ndomain
208        CALL swap_dimensions(ind)
209        CALL swap_geometry(ind)
210        u=f_u(ind)
211        du=f_du(ind)
212        u=du/dumax
213      ENDDO
214      IF (is_mpi_root) PRINT *,"gradiv : it :",it ,": dumax",dumax
215
216    ENDDO 
217    IF (is_mpi_root) PRINT *,"gradiv : dumax",dumax
218    IF (is_mpi_root) PRINT *, 'mean T-cell edge size (km)', 1.45*radius/iim_glo/1000., &
219                              'effective T-cell half-edge size (km)', dumax**(-.5/nitergdiv)/1000
220    IF (is_mpi_root)  PRINT *, 'Max. time step assuming c=340 m/s and Courant number=2.8 :', &
221                               2.8/340.*dumax**(-.5/nitergdiv)
222
223    cgraddiv=dumax**(-1./nitergdiv)
224    IF (is_mpi_root) PRINT *,"cgraddiv : ",cgraddiv
225
226    DO ind=1,ndomain
227      CALL swap_dimensions(ind)
228      CALL swap_geometry(ind)
229      u=f_u(ind)
230
231     DO j=jj_begin,jj_end
232        DO i=ii_begin,ii_end
233          n=(j-1)*iim+i
234          CALL RANDOM_NUMBER(r)
235          u(n+u_right)=r-0.5
236          CALL RANDOM_NUMBER(r)
237          u(n+u_lup)=r-0.5
238          CALL RANDOM_NUMBER(r)
239          u(n+u_ldown)=r-0.5
240       ENDDO
241      ENDDO
242       
243    ENDDO
244
245
246    DO it=1,20
247 
248      dumax=0
249      DO iter=1,nitergrot
250        CALL transfert_request(f_u,req_e1_vect)
251        DO ind=1,ndomain
252          CALL swap_dimensions(ind)
253          CALL swap_geometry(ind)
254          u=f_u(ind)
255          du=f_du(ind)
256          CALL compute_gradrot(u,du,1,1)
257          u=du
258        ENDDO
259      ENDDO
260
261      CALL transfert_request(f_du,req_e1_vect)
262     
263      DO ind=1,ndomain
264        CALL swap_dimensions(ind)
265        CALL swap_geometry(ind)
266        u=f_u(ind)
267        du=f_du(ind)
268       
269        DO j=jj_begin,jj_end
270          DO i=ii_begin,ii_end
271            n=(j-1)*iim+i
272            if (le(n+u_right)>1e-100) dumax=MAX(dumax,ABS(du(n+u_right)))
273            if (le(n+u_lup)>1e-100) dumax=MAX(dumax,ABS(du(n+u_lup)))
274            if (le(n+u_ldown)>1e-100) dumax=MAX(dumax,ABS(du(n+u_ldown)))
275          ENDDO
276        ENDDO
277      ENDDO
278
279      IF (using_mpi) THEN
280        CALL MPI_ALLREDUCE(dumax,dumax1,1,MPI_REAL8,MPI_MAX,comm_icosa,ierr)
281        dumax=dumax1
282      ENDIF
283     
284      DO ind=1,ndomain
285        CALL swap_dimensions(ind)
286        CALL swap_geometry(ind)
287        u=f_u(ind)
288        du=f_du(ind)
289        u=du/dumax
290      ENDDO
291     
292      IF (is_mpi_root) PRINT *,"gradrot : it :",it ,": dumax",dumax
293
294    ENDDO 
295    IF (is_mpi_root) PRINT *,"gradrot : dumax",dumax
296 
297    cgradrot=dumax**(-1./nitergrot) 
298    IF (is_mpi_root) PRINT *,"cgradrot : ",cgradrot
299   
300
301
302    DO ind=1,ndomain
303      CALL swap_dimensions(ind)
304      CALL swap_geometry(ind)
305      theta=f_theta(ind)
306 
307       DO j=jj_begin,jj_end
308          DO i=ii_begin,ii_end
309            n=(j-1)*iim+i
310            CALL RANDOM_NUMBER(r)
311            theta(n)=r-0.5
312         ENDDO
313        ENDDO
314       
315    ENDDO
316
317    DO it=1,20
318
319      dthetamax=0
320      DO iter=1,niterdivgrad
321        CALL transfert_request(f_theta,req_i1)
322        DO ind=1,ndomain
323          CALL swap_dimensions(ind)
324          CALL swap_geometry(ind)
325          theta=f_theta(ind)
326          dtheta=f_dtheta(ind)
327          CALL compute_divgrad(theta,dtheta,1,1)
328          theta=dtheta
329        ENDDO
330      ENDDO
331
332      CALL transfert_request(f_dtheta,req_i1)
333     
334      DO ind=1,ndomain
335        CALL swap_dimensions(ind)
336        CALL swap_geometry(ind)
337        theta=f_theta(ind)
338        dtheta=f_dtheta(ind)
339       
340        DO j=jj_begin,jj_end
341          DO i=ii_begin,ii_end
342            n=(j-1)*iim+i
343            dthetamax=MAX(dthetamax,ABS(dtheta(n)))
344          ENDDO
345        ENDDO
346      ENDDO
347      IF (using_mpi) THEN
348        CALL MPI_ALLREDUCE(dthetamax,dthetamax1,1,MPI_REAL8,MPI_MAX,comm_icosa,ierr)
349        dthetamax=dthetamax1
350      ENDIF
351      IF (is_mpi_root) PRINT *,"divgrad : it :",it ,": dthetamax",dthetamax
352
353      DO ind=1,ndomain
354        CALL swap_dimensions(ind)
355        CALL swap_geometry(ind)
356        theta=f_theta(ind)
357        dtheta=f_dtheta(ind)
358        theta=dtheta/dthetamax
359      ENDDO
360    ENDDO 
361
362    IF (is_mpi_root) PRINT *,"divgrad : divgrad",dthetamax
363 
364    cdivgrad=dthetamax**(-1./niterdivgrad) 
365    IF (is_mpi_root) PRINT *,"cdivgrad : ",cdivgrad
366
367     
368    fact=2
369    DO l=1,llm
370       zz= 1. - preff/presnivs(l)
371       zvert=fact-(fact-1)/(1+zz*zz)
372       tau_graddiv(l) = tau_graddiv(l)/zvert
373       tau_gradrot(l) = tau_gradrot(l)/zvert
374       tau_divgrad(l) = tau_divgrad(l)/zvert
375    ENDDO
376
377    mintau=tau_graddiv(1)
378    DO l=1,llm
379      mintau=MIN(mintau,tau_graddiv(l))
380      mintau=MIN(mintau,tau_gradrot(l))
381      mintau=MIN(mintau,tau_divgrad(l))
382    ENDDO
383       
384    itau_dissip=INT(mintau/dt)
385    itau_dissip=MAX(1,itau_dissip)
386    dtdissip=itau_dissip*dt
387    IF (is_mpi_root) PRINT *,"mintau ",mintau,"itau_dissip",itau_dissip," dtdissip ",dtdissip
388
389  END SUBROUTINE init_dissip 
390 
391 
392  SUBROUTINE dissip(f_ue,f_due,f_ps,f_phis,f_theta_rhodz,f_dtheta_rhodz)
393  USE icosa
394  USE theta2theta_rhodz_mod
395  USE pression_mod
396  USE exner_mod
397  USE geopotential_mod
398  USE trace
399  USE time_mod
400  USE omp_para
401  IMPLICIT NONE
402    TYPE(t_field),POINTER :: f_ue(:)
403    TYPE(t_field),POINTER :: f_due(:)
404    TYPE(t_field),POINTER :: f_ps(:), f_phis(:)
405    TYPE(t_field),POINTER :: f_theta_rhodz(:)
406    TYPE(t_field),POINTER :: f_dtheta_rhodz(:)
407
408    REAL(rstd),POINTER         :: due(:,:)
409    REAL(rstd),POINTER         :: phi(:,:), ue(:,:)
410    REAL(rstd),POINTER         :: due_diss1(:,:)
411    REAL(rstd),POINTER         :: due_diss2(:,:)
412    REAL(rstd),POINTER         :: dtheta_rhodz(:,:)
413    REAL(rstd),POINTER         :: dtheta_rhodz_diss(:,:)
414
415    INTEGER :: ind, shear
416    INTEGER :: l,i,j,n
417
418!$OMP BARRIER
419   
420    CALL trace_start("dissip")
421    CALL gradiv(f_ue,f_due_diss1)
422    CALL gradrot(f_ue,f_due_diss2)
423   
424    CALL divgrad_theta_rhodz(f_ps,f_theta_rhodz,f_dtheta_rhodz_diss)
425   
426! later for openmp   
427!    IF(rayleigh_friction_type>0) THEN
428!       CALL pression(f_ps, f_p)
429!       CALL exner(f_ps, f_p, f_pks, f_pk)
430!       CALL geopotential(f_phis,f_pks,f_pk,f_theta,f_phi)
431!    END IF
432
433    DO ind=1,ndomain
434      CALL swap_dimensions(ind)
435      CALL swap_geometry(ind)
436      due=f_due(ind) 
437      due_diss1=f_due_diss1(ind)
438      due_diss2=f_due_diss2(ind)
439      dtheta_rhodz=f_dtheta_rhodz(ind)
440      dtheta_rhodz_diss=f_dtheta_rhodz_diss(ind)
441
442      DO l=ll_begin,ll_end
443        DO j=jj_begin,jj_end
444          DO i=ii_begin,ii_end
445            n=(j-1)*iim+i
446
447            due(n+u_right,l) = -0.5*( due_diss1(n+u_right,l)/tau_graddiv(l) + due_diss2(n+u_right,l)/tau_gradrot(l))*itau_dissip 
448            due(n+u_lup,l)   = -0.5*( due_diss1(n+u_lup,l)  /tau_graddiv(l) + due_diss2(n+u_lup,l)  /tau_gradrot(l))*itau_dissip
449            due(n+u_ldown,l) = -0.5*( due_diss1(n+u_ldown,l)/tau_graddiv(l) + due_diss2(n+u_ldown,l)/tau_gradrot(l))*itau_dissip
450
451            dtheta_rhodz(n,l) = -0.5*dtheta_rhodz_diss(n,l)/tau_divgrad(l)*itau_dissip
452          ENDDO
453        ENDDO
454      ENDDO
455
456!      dtheta_rhodz=0
457!      due=0
458
459! later for openmp 
460!      IF(rayleigh_friction_type>0) THEN
461!         phi=f_phi(ind)
462!         ue=f_ue(ind)
463!         DO l=1,llm
464!            DO j=jj_begin,jj_end
465!               DO i=ii_begin,ii_end
466!                  n=(j-1)*iim+i
467!                  CALL relax(t_right, u_right)
468!                  CALL relax(t_lup,   u_lup)
469!                  CALL relax(t_ldown, u_ldown)
470!               ENDDO
471!            ENDDO
472!         END DO
473!      END IF
474   END DO
475
476   CALL trace_end("dissip")
477
478!$OMP BARRIER
479
480    CONTAINS
481      SUBROUTINE relax(shift_t, shift_u)
482        USE dcmip_initial_conditions_test_1_2_3
483        REAL(rstd) :: z, ulon,ulat, lon,lat, & ! input to test2_schaer_mountain
484             p,hyam,hybm,w,t,phis,ps,rho,q, &   ! unused input/output to test2_schaer_mountain
485             fz, u3d(3), uref
486        REAL(rstd), PARAMETER :: zh=2e4,ztop=3e4  ! DCMIP values
487        LOGICAL :: hybrid_eta
488        INTEGER :: shift_u, shift_t, zcoords, nn
489        z = (phi(n,l)+phi(n+shift_t,l))/(2.*g)
490        IF(z>zh) THEN  ! relax only in the sponge layer z>zh
491
492           nn = n+shift_u
493           CALL xyz2lonlat(xyz_e(nn,:),lon,lat)
494           zcoords = 1 ; hybrid_eta = .FALSE. ! use z instead of p or hyam/hybm
495           CALL test2_schaer_mountain(lon,lat,p,z,zcoords,hybrid_eta, &
496                hyam,hybm,shear,ulon,ulat,w,t,phis,ps,rho,q)
497!           u3d = ulon*elon_e(nn,:) + ulat*elat_e(nn,:)
498           u3d = ulon*elon_e(nn,:) ! ulat=0
499           uref = sum(u3d*ep_e(nn,:))
500
501           fz = sin((pi/2)*(z-zh)/(ztop-zh))
502           fz = fz*fz/rayleigh_tau
503!           fz = 1/rayleigh_tau
504           due(nn,l) = due(nn,l) - fz*(ue(nn,l)-uref)
505!           due(nn,l) = due(nn,l) - fz*ue(nn,l)
506         END IF
507       END SUBROUTINE relax
508     
509  END SUBROUTINE dissip
510
511  SUBROUTINE gradiv(f_ue,f_due)
512  USE icosa
513  USE trace
514  USE omp_para
515  IMPLICIT NONE
516    TYPE(t_field),POINTER :: f_ue(:)
517    TYPE(t_field),POINTER :: f_due(:)
518    REAL(rstd),POINTER    :: ue(:,:)
519    REAL(rstd),POINTER    :: due(:,:)
520    INTEGER :: ind
521    INTEGER :: it,i,j,l,ij
522       
523    CALL trace_start("gradiv")
524
525    DO ind=1,ndomain
526      CALL swap_dimensions(ind)
527      CALL swap_geometry(ind)
528      ue=f_ue(ind)
529      due=f_due(ind) 
530      DO  l = ll_begin, ll_end
531        DO j=jj_begin,jj_end
532          DO i=ii_begin,ii_end
533            ij=(j-1)*iim+i
534             due(ij+u_right,l)=ue(ij+u_right,l)
535             due(ij+u_lup,l)=ue(ij+u_lup,l)
536             due(ij+u_ldown,l)=ue(ij+u_ldown,l)
537          ENDDO
538        ENDDO
539      ENDDO
540    ENDDO
541
542    DO it=1,nitergdiv
543       
544      CALL send_message(f_due,req_due)
545      CALL wait_message(req_due)
546       
547      DO ind=1,ndomain
548        CALL swap_dimensions(ind)
549        CALL swap_geometry(ind)
550        due=f_due(ind) 
551        CALL compute_gradiv(due,due,ll_begin,ll_end)
552      ENDDO
553    ENDDO
554
555   CALL trace_end("gradiv")
556
557  END SUBROUTINE gradiv
558 
559
560  SUBROUTINE gradrot(f_ue,f_due)
561  USE icosa
562  USE trace
563  USE omp_para
564  IMPLICIT NONE
565    TYPE(t_field),POINTER :: f_ue(:)
566    TYPE(t_field),POINTER :: f_due(:)
567    REAL(rstd),POINTER    :: ue(:,:)
568    REAL(rstd),POINTER    :: due(:,:)
569    INTEGER :: ind
570    INTEGER :: it,i,j,l,ij
571       
572    CALL trace_start("gradrot")
573
574    DO ind=1,ndomain
575      CALL swap_dimensions(ind)
576      CALL swap_geometry(ind)
577      ue=f_ue(ind)
578      due=f_due(ind) 
579      DO  l = ll_begin, ll_end
580        DO j=jj_begin,jj_end
581          DO i=ii_begin,ii_end
582            ij=(j-1)*iim+i
583             due(ij+u_right,l)=ue(ij+u_right,l)
584             due(ij+u_lup,l)=ue(ij+u_lup,l)
585             due(ij+u_ldown,l)=ue(ij+u_ldown,l)
586          ENDDO
587        ENDDO
588      ENDDO
589    ENDDO
590
591    DO it=1,nitergrot
592       
593      CALL send_message(f_due,req_due)
594      CALL wait_message(req_due)
595       
596      DO ind=1,ndomain
597        CALL swap_dimensions(ind)
598        CALL swap_geometry(ind)
599        due=f_due(ind) 
600        CALL compute_gradrot(due,due,ll_begin,ll_end)
601      ENDDO
602
603    ENDDO
604
605    CALL trace_end("gradrot")
606
607  END SUBROUTINE gradrot
608 
609  SUBROUTINE divgrad(f_theta,f_dtheta)
610  USE icosa
611  USE trace
612  USE omp_para
613  IMPLICIT NONE
614    TYPE(t_field),POINTER :: f_theta(:)
615    TYPE(t_field),POINTER :: f_dtheta(:)
616    REAL(rstd),POINTER    :: theta(:,:)
617    REAL(rstd),POINTER    :: dtheta(:,:)
618    INTEGER :: ind
619    INTEGER :: it
620
621    CALL trace_start("divgrad")
622       
623    DO ind=1,ndomain
624      CALL swap_dimensions(ind)
625      CALL swap_geometry(ind)
626      theta=f_theta(ind)
627      dtheta=f_dtheta(ind) 
628      dtheta=theta
629    ENDDO
630
631    DO it=1,niterdivgrad
632       
633      CALL transfert_request(f_dtheta,req_i1)
634       
635      DO ind=1,ndomain
636        CALL swap_dimensions(ind)
637        CALL swap_geometry(ind)
638        dtheta=f_dtheta(ind) 
639        CALL compute_divgrad(dtheta,dtheta,ll_begin,ll_end)
640      ENDDO
641
642    ENDDO
643
644    CALL trace_end("divgrad")
645
646  END SUBROUTINE divgrad
647   
648  SUBROUTINE divgrad_theta_rhodz(f_ps,f_theta_rhodz,f_dtheta_rhodz)
649  USE icosa
650  USE trace
651  USE omp_para
652  USE disvert_mod
653  IMPLICIT NONE
654    TYPE(t_field),POINTER :: f_ps(:)
655    TYPE(t_field),POINTER :: f_theta_rhodz(:)
656    TYPE(t_field),POINTER :: f_dtheta_rhodz(:)
657   
658    REAL(rstd),POINTER :: ps(:)
659    REAL(rstd),POINTER :: theta_rhodz(:,:)
660    REAL(rstd),POINTER :: dtheta_rhodz(:,:)
661
662    INTEGER :: ind
663    INTEGER :: it,i,j,l,ij
664
665    CALL trace_start("divgrad")
666       
667    DO ind=1,ndomain
668      CALL swap_dimensions(ind)
669      CALL swap_geometry(ind)
670      ps=f_ps(ind)
671      theta_rhodz=f_theta_rhodz(ind)
672      dtheta_rhodz=f_dtheta_rhodz(ind) 
673      DO  l = ll_begin, ll_end
674        DO j=jj_begin,jj_end
675          DO i=ii_begin,ii_end
676            ij=(j-1)*iim+i
677            dtheta_rhodz(ij,l) = theta_rhodz(ij,l) / ( ap(l)-ap(l+1) + (bp(l)-bp(l+1))*ps(ij))/g
678          ENDDO
679        ENDDO
680      ENDDO
681    ENDDO
682
683    DO it=1,niterdivgrad
684       
685      CALL send_message(f_dtheta_rhodz,req_dtheta)
686      CALL wait_message(req_dtheta)
687      DO ind=1,ndomain
688        CALL swap_dimensions(ind)
689        CALL swap_geometry(ind)
690        dtheta_rhodz=f_dtheta_rhodz(ind) 
691        CALL compute_divgrad(dtheta_rhodz,dtheta_rhodz,ll_begin,ll_end)
692      ENDDO
693
694    ENDDO
695
696    DO ind=1,ndomain
697      CALL swap_dimensions(ind)
698      CALL swap_geometry(ind)
699      dtheta_rhodz=f_dtheta_rhodz(ind) 
700   
701      DO  l = ll_begin, ll_end
702        DO j=jj_begin,jj_end
703          DO i=ii_begin,ii_end
704            ij=(j-1)*iim+i
705            dtheta_rhodz(ij,l) = dtheta_rhodz(ij,l) * ( ap(l)-ap(l+1) + (bp(l)-bp(l+1))*ps(ij))/g
706          ENDDO
707        ENDDO
708      ENDDO
709    ENDDO
710
711
712    CALL trace_end("divgrad")
713
714  END SUBROUTINE divgrad_theta_rhodz 
715 
716  SUBROUTINE compute_gradiv(ue,gradivu_e,llb,lle)
717  USE icosa
718  IMPLICIT NONE
719    INTEGER,INTENT(IN)     :: llb
720    INTEGER,INTENT(IN)     :: lle
721    REAL(rstd),INTENT(IN)  :: ue(iim*3*jjm,llm)
722    REAL(rstd),INTENT(OUT) :: gradivu_e(iim*3*jjm,llm)
723    REAL(rstd) :: divu_i(iim*jjm,llb:lle)
724   
725    INTEGER :: i,j,n,l
726
727    DO l=llb,lle
728      DO j=jj_begin,jj_end
729        DO i=ii_begin,ii_end
730          n=(j-1)*iim+i
731          divu_i(n,l)=1./Ai(n)*(ne(n,right)*ue(n+u_right,l)*le(n+u_right)  +  &
732                             ne(n,rup)*ue(n+u_rup,l)*le(n+u_rup)        +  & 
733                             ne(n,lup)*ue(n+u_lup,l)*le(n+u_lup)        +  & 
734                             ne(n,left)*ue(n+u_left,l)*le(n+u_left)     +  & 
735                             ne(n,ldown)*ue(n+u_ldown,l)*le(n+u_ldown)  +  & 
736                             ne(n,rdown)*ue(n+u_rdown,l)*le(n+u_rdown))
737        ENDDO
738      ENDDO
739    ENDDO
740   
741    DO l=llb,lle
742      DO j=jj_begin,jj_end
743        DO i=ii_begin,ii_end
744 
745          n=(j-1)*iim+i
746 
747          gradivu_e(n+u_right,l)=-1/de(n+u_right)*(ne(n,right)*divu_i(n,l)+ ne(n+t_right,left)*divu_i(n+t_right,l) )       
748
749          gradivu_e(n+u_lup,l)=-1/de(n+u_lup)*(ne(n,lup)*divu_i(n,l)+ ne(n+t_lup,rdown)*divu_i(n+t_lup,l))       
750   
751          gradivu_e(n+u_ldown,l)=-1/de(n+u_ldown)*(ne(n,ldown)*divu_i(n,l)+ne(n+t_ldown,rup)*divu_i(n+t_ldown,l) )
752
753        ENDDO
754      ENDDO
755    ENDDO
756
757    DO l=llb,lle
758      DO j=jj_begin,jj_end
759        DO i=ii_begin,ii_end
760          n=(j-1)*iim+i
761          gradivu_e(n+u_right,l)=-gradivu_e(n+u_right,l)*cgraddiv
762          gradivu_e(n+u_lup,l)=-gradivu_e(n+u_lup,l)*cgraddiv
763          gradivu_e(n+u_ldown,l)=-gradivu_e(n+u_ldown,l)*cgraddiv
764        ENDDO
765      ENDDO
766    ENDDO
767
768   
769  END SUBROUTINE compute_gradiv
770 
771  SUBROUTINE compute_divgrad(theta,divgrad_i,llb,lle)
772  USE icosa
773  IMPLICIT NONE
774    INTEGER,INTENT(IN)     :: llb
775    INTEGER,INTENT(IN)     :: lle
776    REAL(rstd),INTENT(IN)  :: theta(iim*jjm,llm)
777    REAL(rstd),INTENT(OUT) :: divgrad_i(iim*jjm,llm)
778    REAL(rstd)  :: grad_e(3*iim*jjm,llb:lle)
779
780    INTEGER :: i,j,n,l
781
782       
783    DO l=llb,lle
784      DO j=jj_begin-1,jj_end+1
785        DO i=ii_begin-1,ii_end+1
786   
787          n=(j-1)*iim+i
788 
789          grad_e(n+u_right,l)=-1/de(n+u_right)*(ne(n,right)*theta(n,l)+ ne(n+t_right,left)*theta(n+t_right,l) )       
790
791          grad_e(n+u_lup,l)=-1/de(n+u_lup)*(ne(n,lup)*theta(n,l)+ ne(n+t_lup,rdown)*theta(n+t_lup,l ))       
792   
793          grad_e(n+u_ldown,l)=-1/de(n+u_ldown)*(ne(n,ldown)*theta(n,l)+ne(n+t_ldown,rup)*theta(n+t_ldown,l) )
794
795        ENDDO
796      ENDDO
797    ENDDO
798   
799   
800    DO l=llb,lle
801      DO j=jj_begin,jj_end
802        DO i=ii_begin,ii_end
803          n=(j-1)*iim+i
804          divgrad_i(n,l)=1./Ai(n)*(ne(n,right)*grad_e(n+u_right,l)*le(n+u_right)  +  &
805                             ne(n,rup)*grad_e(n+u_rup,l)*le(n+u_rup)              +  & 
806                             ne(n,lup)*grad_e(n+u_lup,l)*le(n+u_lup)              +  & 
807                             ne(n,left)*grad_e(n+u_left,l)*le(n+u_left)           +  & 
808                             ne(n,ldown)*grad_e(n+u_ldown,l)*le(n+u_ldown)        +  & 
809                             ne(n,rdown)*grad_e(n+u_rdown,l)*le(n+u_rdown))
810        ENDDO
811      ENDDO
812    ENDDO
813   
814    DO l=llb,lle
815      DO j=jj_begin,jj_end
816        DO i=ii_begin,ii_end
817          n=(j-1)*iim+i
818          divgrad_i(n,l)=-divgrad_i(n,l)*cdivgrad
819        ENDDO
820      ENDDO
821    ENDDO
822
823  END SUBROUTINE compute_divgrad
824
825   
826  SUBROUTINE compute_gradrot(ue,gradrot_e,llb,lle)
827  USE icosa
828  IMPLICIT NONE
829    INTEGER,INTENT(IN)     :: llb
830    INTEGER,INTENT(IN)     :: lle
831    REAL(rstd),INTENT(IN)  :: ue(iim*3*jjm,llm)
832    REAL(rstd),INTENT(OUT) :: gradrot_e(iim*3*jjm,llm)
833    REAL(rstd) :: rot_v(2*iim*jjm,llb:lle)
834
835    INTEGER :: i,j,n,l
836     
837    DO l=llb,lle
838      DO j=jj_begin-1,jj_end+1
839        DO i=ii_begin-1,ii_end+1
840          n=(j-1)*iim+i
841       
842          rot_v(n+z_up,l)= 1./Av(n+z_up)*(  ne(n,rup)*ue(n+u_rup,l)*de(n+u_rup)                     &
843                                + ne(n+t_rup,left)*ue(n+t_rup+u_left,l)*de(n+t_rup+u_left)          &
844                                - ne(n,lup)*ue(n+u_lup,l)*de(n+u_lup) ) 
845                             
846          rot_v(n+z_down,l) = 1./Av(n+z_down)*( ne(n,ldown)*ue(n+u_ldown,l)*de(n+u_ldown)                 &
847                                     + ne(n+t_ldown,right)*ue(n+t_ldown+u_right,l)*de(n+t_ldown+u_right)  &
848                                     - ne(n,rdown)*ue(n+u_rdown,l)*de(n+u_rdown) )
849 
850        ENDDO
851      ENDDO
852    ENDDO                             
853 
854    DO l=llb,lle
855      DO j=jj_begin,jj_end
856        DO i=ii_begin,ii_end
857          n=(j-1)*iim+i
858 
859          gradrot_e(n+u_right,l)=1/le(n+u_right)*ne(n,right)*(rot_v(n+z_rdown,l)-rot_v(n+z_rup,l)) 
860         
861          gradrot_e(n+u_lup,l)=1/le(n+u_lup)*ne(n,lup)*(rot_v(n+z_up,l)-rot_v(n+z_lup,l)) 
862       
863          gradrot_e(n+u_ldown,l)=1/le(n+u_ldown)*ne(n,ldown)*(rot_v(n+z_ldown,l)-rot_v(n+z_down,l))
864       
865        ENDDO
866      ENDDO
867    ENDDO
868
869    DO l=llb,lle
870      DO j=jj_begin,jj_end
871        DO i=ii_begin,ii_end
872          n=(j-1)*iim+i
873   
874          gradrot_e(n+u_right,l)=-gradrot_e(n+u_right,l)*cgradrot       
875          gradrot_e(n+u_lup,l)=-gradrot_e(n+u_lup,l)*cgradrot
876          gradrot_e(n+u_ldown,l)=-gradrot_e(n+u_ldown,l)*cgradrot
877       
878        ENDDO
879      ENDDO
880    ENDDO 
881   
882  END SUBROUTINE compute_gradrot
883
884
885END MODULE dissip_gcm_mod
886       
Note: See TracBrowser for help on using the repository browser.