source: codes/icosagcm/branches/SATURN_DYNAMICO/ICOSAGCM/src/physics_lmdz_generic.f90 @ 246

Last change on this file since 246 was 246, checked in by ymipsl, 10 years ago

One call for initialize physics from dynamico.

YM

File size: 16.4 KB
Line 
1MODULE physics_lmdz_generic_mod
2  USE field_mod
3  USE transfert_mod
4 
5  INTEGER,SAVE :: nbp_phys
6  TYPE(t_message) :: req_u
7
8  TYPE(t_field),POINTER :: f_p(:) 
9  TYPE(t_field),POINTER :: f_pks(:) 
10  TYPE(t_field),POINTER :: f_pk(:) 
11  TYPE(t_field),POINTER :: f_p_layer(:)   
12  TYPE(t_field),POINTER :: f_theta(:)   
13  TYPE(t_field),POINTER :: f_phi(:)   
14  TYPE(t_field),POINTER :: f_Temp(:)   
15  TYPE(t_field),POINTER :: f_ulon(:)   
16  TYPE(t_field),POINTER :: f_ulat(:)   
17  TYPE(t_field),POINTER :: f_dulon(:)
18  TYPE(t_field),POINTER :: f_dulat(:)
19  TYPE(t_field),POINTER :: f_dTemp(:)
20  TYPE(t_field),POINTER :: f_dq(:)
21  TYPE(t_field),POINTER :: f_dps(:)
22  TYPE(t_field),POINTER :: f_duc(:)
23
24  INTEGER :: start_clock
25  INTEGER :: stop_clock
26  INTEGER :: count_clock=0
27 
28  REAL :: start_day
29  REAL :: day_length
30  REAL :: year_length
31 
32 INTEGER,ALLOCATABLE,SAVE :: domain_offset(:)
33
34 
35
36CONTAINS
37
38  SUBROUTINE init_physics
39  USE icosa
40  USE domain_mod
41  USE dimensions
42  USE mpi_mod
43  USE mpipara
44  USE disvert_mod
45 
46  IMPLICIT NONE
47  INTEGER  :: distrib(0:mpi_size-1)
48  INTEGER  :: ind,i,j,ij,pos
49
50  REAL(rstd),ALLOCATABLE :: latfi(:)
51  REAL(rstd),ALLOCATABLE :: lonfi(:)
52  REAL(rstd),ALLOCATABLE :: airefi(:)
53
54    start_day=0
55    day_length=86400
56    year_length=86400*365.25
57   
58    CALL getin('start_day',start_day)
59    CALL getin('day_length',day_length)
60    CALL getin('year_length',year_length)
61
62!$OMP PARALLEL
63    CALL allocate_field(f_p,field_t,type_real,llm+1)
64    CALL allocate_field(f_pks,field_t,type_real)
65    CALL allocate_field(f_pk,field_t,type_real,llm)
66    CALL allocate_field(f_p_layer,field_t,type_real,llm)
67    CALL allocate_field(f_theta,field_t,type_real,llm)
68    CALL allocate_field(f_phi,field_t,type_real,llm)
69    CALL allocate_field(f_Temp,field_t,type_real,llm)
70    CALL allocate_field(f_ulon,field_t,type_real,llm)
71    CALL allocate_field(f_ulat,field_t,type_real,llm)
72    CALL allocate_field(f_dulon,field_t,type_real,llm)
73    CALL allocate_field(f_dulat,field_t,type_real,llm)
74    CALL allocate_field(f_dTemp,field_t,type_real,llm)
75    CALL allocate_field(f_dq,field_t,type_real,llm,nqtot)
76    CALL allocate_field(f_dps,field_t,type_real)
77    CALL allocate_field(f_duc,field_t,type_real,3,llm)   
78!$OMP END PARALLEL   
79
80    ALLOCATE(domain_offset(ndomain))
81    nbp_phys=0
82    domain_offset(1)=0
83    DO ind=1,ndomain
84      CALL swap_dimensions(ind)
85      IF (ind<ndomain) THEN
86        domain_offset(ind+1)=domain_offset(ind)+ii_nb*jj_nb
87      ENDIF
88      nbp_phys=nbp_phys+ii_nb*jj_nb
89    ENDDO
90   
91    CALL MPI_ALLGATHER(nbp_phys,1,MPI_INTEGER,distrib,1,MPI_INTEGER,comm_icosa,ierr)
92   
93    ALLOCATE(latfi(nbp_phys))
94    ALLOCATE(lonfi(nbp_phys))
95    ALLOCATE(airefi(nbp_phys))
96   
97    pos=0
98    DO ind=1,ndomain
99      CALL swap_dimensions(ind)
100      CALL swap_geometry(ind)
101      DO j=jj_begin,jj_end
102        DO i=ii_begin,ii_end
103          ij=(j-1)*iim+i
104          pos=pos+1
105          CALL xyz2lonlat(xyz_i(ij,:),lonfi(pos),latfi(pos))
106          airefi(pos)=Ai(ij) 
107        ENDDO
108      ENDDO
109    ENDDO
110
111
112    CALL initialize_unstructured_physics(nbp_phys,llm, comm_icosa, mpi_size,distrib,            &
113                                         day_length,start_day,itau_physics*dt,                  &
114                                         latfi,lonfi,airefi,radius,g, gas_constant/mu, cpp,     &
115                                         preff, ap, bp)                                         
116   
117   
118   
119!    CALL init_phys_lmdz(128,97,llm, comm_icosa, mpi_size, distrib)
120   
121    CALL init_gcm_lmdz(nbp_phys,mpi_size,distrib,latfi,lonfi,airefi)
122!    CALL init_phys_lmdz(128,97,llm, comm_icosa, mpi_size, distrib)
123!    CALL iniphysiq(llm,day_length,start_day,itau_physics*dt,     &
124!                   latfi,lonfi,airefi,radius,g, gas_constant/mu, cpp,     &
125!                   ap, bp)
126   
127   
128   
129  END SUBROUTINE init_physics
130 
131  SUBROUTINE physics(it,f_phis, f_ps, f_theta_rhodz, f_u, f_wflux, f_q)
132  USE ICOSA
133  USE time_mod
134  USE disvert_mod
135  USE transfert_mod
136  USE mpipara
137  IMPLICIT NONE
138    INTEGER,INTENT(IN)    :: it
139    TYPE(t_field),POINTER :: f_phis(:)
140    TYPE(t_field),POINTER :: f_ps(:)
141    TYPE(t_field),POINTER :: f_theta_rhodz(:)
142    TYPE(t_field),POINTER :: f_u(:)
143    TYPE(t_field),POINTER :: f_wflux(:)
144    TYPE(t_field),POINTER :: f_q(:)
145 
146    REAL(rstd),POINTER :: phis(:)
147    REAL(rstd),POINTER :: ps(:)
148    REAL(rstd),POINTER :: theta_rhodz(:,:)
149    REAL(rstd),POINTER :: u(:,:)
150    REAL(rstd),POINTER :: wflux(:,:)
151    REAL(rstd),POINTER :: q(:,:,:)
152    REAL(rstd),POINTER :: p(:,:)
153    REAL(rstd),POINTER :: pks(:)
154    REAL(rstd),POINTER :: pk(:,:)
155    REAL(rstd),POINTER :: p_layer(:,:)
156    REAL(rstd),POINTER :: theta(:,:)
157    REAL(rstd),POINTER :: phi(:,:)
158    REAL(rstd),POINTER :: Temp(:,:)
159    REAL(rstd),POINTER :: ulon(:,:)
160    REAL(rstd),POINTER :: ulat(:,:)
161    REAL(rstd),POINTER :: dulon(:,:)
162    REAL(rstd),POINTER :: dulat(:,:)
163    REAL(rstd),POINTER :: dTemp(:,:)
164    REAL(rstd),POINTER :: dq(:,:,:)
165    REAL(rstd),POINTER :: dps(:)
166    REAL(rstd),POINTER :: duc(:,:,:)
167
168
169    INTEGER :: ind
170   
171    REAL(rstd),ALLOCATABLE,SAVE :: ps_phy(:)
172    REAL(rstd),ALLOCATABLE,SAVE :: p_phy(:,:)
173    REAL(rstd),ALLOCATABLE,SAVE :: p_layer_phy(:,:)
174    REAL(rstd),ALLOCATABLE,SAVE :: Temp_phy(:,:)
175    REAL(rstd),ALLOCATABLE,SAVE :: phis_phy(:)
176    REAL(rstd),ALLOCATABLE,SAVE :: phi_phy(:,:)
177    REAL(rstd),ALLOCATABLE,SAVE :: ulon_phy(:,:)
178    REAL(rstd),ALLOCATABLE,SAVE :: ulat_phy(:,:)
179    REAL(rstd),ALLOCATABLE,SAVE :: q_phy(:,:,:)
180    REAL(rstd),ALLOCATABLE,SAVE :: wflux_phy(:,:)
181    REAL(rstd),ALLOCATABLE,SAVE :: dulon_phy(:,:)
182    REAL(rstd),ALLOCATABLE,SAVE :: dulat_phy(:,:)
183    REAL(rstd),ALLOCATABLE,SAVE :: dTemp_phy(:,:)
184    REAL(rstd),ALLOCATABLE,SAVE :: dq_phy(:,:,:)
185    REAL(rstd),ALLOCATABLE,SAVE :: dps_phy(:)
186    REAL(rstd) :: dtphy 
187    REAL(rstd) :: time
188    REAL(rstd) :: day
189    REAL(rstd) :: real_time
190    INTEGER    :: offset
191    LOGICAL :: lafin=.FALSE.
192    LOGICAL,SAVE :: first=.TRUE.
193!$OMP THREADPRIVATE(first)
194
195   
196    IF (first) THEN
197      first=.FALSE.
198      CALL init_message(f_u,req_e1_vect,req_u)
199!$OMP MASTER
200      ALLOCATE(ps_phy(nbp_phys))
201      ALLOCATE(p_phy(nbp_phys,llm+1))
202      ALLOCATE(p_layer_phy(nbp_phys,llm))
203      ALLOCATE(Temp_phy(nbp_phys,llm))
204      ALLOCATE(phis_phy(nbp_phys))
205      ALLOCATE(phi_phy(nbp_phys,llm))
206      ALLOCATE(ulon_phy(nbp_phys,llm))
207      ALLOCATE(ulat_phy(nbp_phys,llm))
208      ALLOCATE(q_phy(nbp_phys,llm,nqtot))
209      ALLOCATE(wflux_phy(nbp_phys,llm))
210      ALLOCATE(dulon_phy(nbp_phys,llm))
211      ALLOCATE(dulat_phy(nbp_phys,llm))
212      ALLOCATE(dTemp_phy(nbp_phys,llm))
213      ALLOCATE(dq_phy(nbp_phys,llm,nqtot))
214      ALLOCATE(dps_phy(nbp_phys))
215     
216!$OMP END MASTER
217!$OMP BARRIER
218    ENDIF
219
220!$OMP MASTER       
221!    CALL update_calendar(it)
222!$OMP END MASTER
223!$OMP BARRIER
224    dtphy=itau_physics*dt
225    real_time=start_day*day_length+it*dt
226    day  = INT( MODULO(real_time,year_length) / day_length) 
227    time = MODULO(real_time,day_length) / day_length
228
229!$OMP MASTER   
230    IF (is_mpi_root) PRINT*,"Enterring in physic : day", day, "  time : ",time 
231!$OMP END MASTER   
232   
233   
234   
235   
236    CALL transfert_message(f_u,req_u)
237   
238    DO ind=1,ndomain
239      CALL swap_dimensions(ind)
240      IF (assigned_domain(ind)) THEN
241        offset=domain_offset(ind)
242        CALL swap_geometry(ind)
243     
244        phis=f_phis(ind)
245        ps=f_ps(ind)
246        theta_rhodz=f_theta_rhodz(ind)
247        u=f_u(ind)
248        q=f_q(ind)
249        wflux=f_wflux(ind)
250        p=f_p(ind)
251        pks=f_pks(ind)
252        pk=f_pk(ind)
253        p_layer=f_p_layer(ind)
254        theta=f_theta(ind)
255        phi=f_phi(ind)
256        Temp=f_Temp(ind)
257        ulon=f_ulon(ind)
258        ulat=f_ulat(ind)
259           
260        CALL grid_icosa_to_physics
261
262      ENDIF
263    ENDDO
264
265!$OMP BARRIER
266!$OMP MASTER
267    CALL SYSTEM_CLOCK(start_clock)
268!$OMP END MASTER
269    CALL calfis_icosa(dtphy, lafin, day, time, presnivs,      &
270                      p_phy, p_layer_phy, phi_phy, phis_phy, ulon_phy, ulat_phy, Temp_phy, q_phy, wflux_phy, &
271                      dulon_phy, dulat_phy, dTemp_phy, dq_phy, dps_phy  )
272
273
274!$OMP MASTER
275    CALL SYSTEM_CLOCK(stop_clock)
276    count_clock=count_clock+stop_clock-start_clock
277!$OMP END MASTER
278
279!$OMP BARRIER                       
280
281    DO ind=1,ndomain
282      CALL swap_dimensions(ind)
283      IF (assigned_domain(ind)) THEN
284        CALL swap_geometry(ind)
285        offset=domain_offset(ind)
286     
287        theta_rhodz=f_theta_rhodz(ind)
288        u=f_u(ind)
289        q=f_q(ind)
290        ps=f_ps(ind)
291        dulon=f_dulon(ind)
292        dulat=f_dulat(ind)
293        Temp=f_temp(ind)
294        dTemp=f_dTemp(ind)
295        dq=f_dq(ind)
296        dps=f_dps(ind)
297        duc=f_duc(ind)
298        p=f_p(ind)
299        pks=f_pks(ind)
300        pk=f_pk(ind)
301     
302        CALL grid_physics_to_icosa
303      ENDIF
304    ENDDO
305
306!$OMP BARRIER
307     
308  CONTAINS
309   
310    SUBROUTINE grid_physics_to_icosa
311    USE theta2theta_rhodz_mod
312    USE omp_para
313    IMPLICIT NONE
314      INTEGER :: i,j,ij,l,iq
315         
316      DO j=jj_begin,jj_end
317        DO i=ii_begin,ii_end
318          ij=(j-1)*iim+i
319          offset=offset+1
320         
321          dulon(ij,ll_begin:ll_end)=dulon_phy(offset,ll_begin:ll_end)
322          dulat(ij,ll_begin:ll_end)=dulat_phy(offset,ll_begin:ll_end)
323          dTemp(ij,ll_begin:ll_end)=dTemp_phy(offset,ll_begin:ll_end)
324          Temp(ij,ll_begin:ll_end) = Temp_phy(offset,ll_begin:ll_end)
325          dq(ij,ll_begin:ll_end,:)=dq_phy(offset,ll_begin:ll_end,:)
326          if (omp_first) dps(ij)=dps_phy(offset)
327        ENDDO
328      ENDDO
329   
330      DO l=ll_begin,ll_end
331        DO j=jj_begin,jj_end
332          DO i=ii_begin,ii_end
333            ij=(j-1)*iim+i
334            duc(ij,:,l)=dulon(ij,l)*elon_i(ij,:)+dulat(ij,l)*elat_i(ij,:)
335          ENDDO
336        ENDDO
337      ENDDO
338
339      DO l=ll_begin,ll_end
340        DO j=jj_begin,jj_end
341          DO i=ii_begin,ii_end
342            ij=(j-1)*iim+i
343            u(ij+u_right,l) = u(ij+u_right,l) + dtphy * sum( 0.5*(duc(ij,:,l) + duc(ij+t_right,:,l))*ep_e(ij+u_right,:) )
344            u(ij+u_lup,l) = u(ij+u_lup,l) + dtphy * sum( 0.5*(duc(ij,:,l) + duc(ij+t_lup,:,l))*ep_e(ij+u_lup,:) )
345            u(ij+u_ldown,l) = u(ij+u_ldown,l) + dtphy*sum( 0.5*(duc(ij,:,l) + duc(ij+t_ldown,:,l))*ep_e(ij+u_ldown,:) )
346          ENDDO
347        ENDDO
348      ENDDO         
349
350      DO l=ll_begin,ll_end
351        DO j=jj_begin,jj_end
352          DO i=ii_begin,ii_end
353            ij=(j-1)*iim+i
354            Temp(ij,l)=Temp(ij,l)+ dtphy * dTemp(ij,l)
355          ENDDO
356        ENDDO
357      ENDDO         
358     
359      DO iq=1,nqtot
360        DO l=ll_begin,ll_end
361          DO j=jj_begin,jj_end
362            DO i=ii_begin,ii_end
363              ij=(j-1)*iim+i
364              q(ij,l,iq)=q(ij,l,iq)+ dtphy * dq(ij,l,iq)
365            ENDDO
366          ENDDO
367        ENDDO 
368      ENDDO
369
370!$OMP BARRIER
371     
372     IF (omp_first) THEN
373       DO j=jj_begin,jj_end
374         DO i=ii_begin,ii_end
375           ij=(j-1)*iim+i
376           ps(ij)=ps(ij)+ dtphy * dps(ij)
377          ENDDO
378       ENDDO
379     ENDIF
380
381!     CALL compute_temperature2theta_rhodz(ps,Temp,theta_rhodz,0)
382
383! compute pression
384!$OMP BARRIER
385      DO    l    = ll_begin,ll_endp1
386        DO j=jj_begin,jj_end
387          DO i=ii_begin,ii_end
388            ij=(j-1)*iim+i
389            p(ij,l) = ap(l) + bp(l) * ps(ij)
390          ENDDO
391        ENDDO
392      ENDDO
393
394!$OMP BARRIER
395
396! compute exner
397       
398       IF (omp_first) THEN
399         DO j=jj_begin,jj_end
400            DO i=ii_begin,ii_end
401               ij=(j-1)*iim+i
402               pks(ij) = cpp * ( ps(ij)/preff ) ** kappa
403            ENDDO
404         ENDDO
405       ENDIF
406
407       ! 3D : pk
408       DO l = ll_begin,ll_end
409          DO j=jj_begin,jj_end
410             DO i=ii_begin,ii_end
411                ij=(j-1)*iim+i
412                pk(ij,l) = cpp * ((.5/preff)*(p(ij,l)+p(ij,l+1))) ** kappa
413             ENDDO
414          ENDDO
415       ENDDO
416
417!$OMP BARRIER
418
419!   compute theta, temperature and pression at layer
420    DO    l    = ll_begin, ll_end
421      DO j=jj_begin,jj_end
422        DO i=ii_begin,ii_end
423          ij=(j-1)*iim+i
424          theta_rhodz(ij,l) = temp(ij,l) * ((p(ij,l)-p(ij,l+1))/g) / (pk(ij,l) / cpp )
425        ENDDO
426      ENDDO
427    ENDDO
428   
429    END SUBROUTINE grid_physics_to_icosa
430   
431   
432    SUBROUTINE grid_icosa_to_physics
433    USE pression_mod
434    USE exner_mod
435    USE theta2theta_rhodz_mod
436    USE geopotential_mod
437    USE wind_mod
438    USE omp_para
439    IMPLICIT NONE
440   
441    REAL(rstd) :: uc(3)
442    INTEGER :: i,j,ij,l
443   
444
445! compute pression
446
447      DO    l    = ll_begin,ll_endp1
448        DO j=jj_begin,jj_end
449          DO i=ii_begin,ii_end
450            ij=(j-1)*iim+i
451            p(ij,l) = ap(l) + bp(l) * ps(ij)
452          ENDDO
453        ENDDO
454      ENDDO
455
456!$OMP BARRIER
457
458! compute exner
459       
460       IF (omp_first) THEN
461         DO j=jj_begin,jj_end
462            DO i=ii_begin,ii_end
463               ij=(j-1)*iim+i
464               pks(ij) = cpp * ( ps(ij)/preff ) ** kappa
465            ENDDO
466         ENDDO
467       ENDIF
468
469       ! 3D : pk
470       DO l = ll_begin,ll_end
471          DO j=jj_begin,jj_end
472             DO i=ii_begin,ii_end
473                ij=(j-1)*iim+i
474                pk(ij,l) = cpp * ((.5/preff)*(p(ij,l)+p(ij,l+1))) ** kappa
475             ENDDO
476          ENDDO
477       ENDDO
478
479!$OMP BARRIER
480
481!   compute theta, temperature and pression at layer
482    DO    l    = ll_begin, ll_end
483      DO j=jj_begin,jj_end
484        DO i=ii_begin,ii_end
485          ij=(j-1)*iim+i
486          theta(ij,l) = theta_rhodz(ij,l) / ((p(ij,l)-p(ij,l+1))/g)
487          Temp(ij,l) = theta(ij,l) * pk(ij,l) / cpp
488          p_layer(ij,l)=preff*(pk(ij,l)/cpp)**(1./kappa) 
489        ENDDO
490      ENDDO
491    ENDDO
492
493
494!!! Compute geopotential
495       
496  ! for first layer
497  IF (omp_first) THEN
498    DO j=jj_begin,jj_end
499      DO i=ii_begin,ii_end
500        ij=(j-1)*iim+i
501        phi( ij,1 ) = phis( ij ) + theta(ij,1) * ( pks(ij) - pk(ij,1) )
502      ENDDO
503    ENDDO
504  ENDIF
505!!-> implicit flush on phi(:,1)
506         
507  ! for other layers
508   DO l = ll_beginp1, ll_end
509     DO j=jj_begin,jj_end
510       DO i=ii_begin,ii_end
511         ij=(j-1)*iim+i
512         phi(ij,l) =  0.5 * ( theta(ij,l)  + theta(ij,l-1) )  & 
513                          * (  pk(ij,l-1) -  pk(ij,l)    )
514       ENDDO
515     ENDDO
516   ENDDO       
517
518!$OMP BARRIER
519
520   DO l = 2, llm
521     DO j=jj_begin,jj_end
522! ---> Bug compilo intel ici en openmp
523! ---> Couper la boucle
524       if (j==jj_end+1) PRINT*,"this message must not be printed"
525       DO i=ii_begin,ii_end
526         ij=(j-1)*iim+i
527         phi(ij,l) = phi(ij,l)+ phi(ij,l-1)
528       ENDDO
529     ENDDO
530   ENDDO
531! --> IMPLICIT FLUSH on phi
532
533
534
535! compute wind centered lon lat compound
536    DO l=ll_begin,ll_end
537      DO j=jj_begin,jj_end
538        DO i=ii_begin,ii_end
539          ij=(j-1)*iim+i
540          uc(:)=1/Ai(ij)*                                                                                                &
541                        ( ne(ij,right)*u(ij+u_right,l)*le(ij+u_right)*((xyz_v(ij+z_rdown,:)+xyz_v(ij+z_rup,:))/2-centroid(ij,:))  &
542                         + ne(ij,rup)*u(ij+u_rup,l)*le(ij+u_rup)*((xyz_v(ij+z_rup,:)+xyz_v(ij+z_up,:))/2-centroid(ij,:))          &
543                         + ne(ij,lup)*u(ij+u_lup,l)*le(ij+u_lup)*((xyz_v(ij+z_up,:)+xyz_v(ij+z_lup,:))/2-centroid(ij,:))          &
544                         + ne(ij,left)*u(ij+u_left,l)*le(ij+u_left)*((xyz_v(ij+z_lup,:)+xyz_v(ij+z_ldown,:))/2-centroid(ij,:))    &
545                         + ne(ij,ldown)*u(ij+u_ldown,l)*le(ij+u_ldown)*((xyz_v(ij+z_ldown,:)+xyz_v(ij+z_down,:))/2-centroid(ij,:))&
546                         + ne(ij,rdown)*u(ij+u_rdown,l)*le(ij+u_rdown)*((xyz_v(ij+z_down,:)+xyz_v(ij+z_rdown,:))/2-centroid(ij,:)))
547          ulon(ij,l)=sum(uc(:)*elon_i(ij,:))
548          ulat(ij,l)=sum(uc(:)*elat_i(ij,:)) 
549        ENDDO
550      ENDDO
551    ENDDO
552
553!$OMP BARRIER
554
555
556     
557      DO j=jj_begin,jj_end
558        DO i=ii_begin,ii_end
559          ij=(j-1)*iim+i
560          offset=offset+1
561
562          IF (omp_first) ps_phy(offset) = ps(ij)
563          p_phy(offset,ll_begin:ll_endp1) = p(ij,ll_begin:ll_endp1)
564          p_layer_phy(offset,ll_begin:ll_end) = p_layer(ij,ll_begin:ll_end)
565          Temp_phy(offset,ll_begin:ll_end) = Temp(ij,ll_begin:ll_end)
566          IF (omp_first) phis_phy(offset) = phis(ij)
567          phi_phy(offset,ll_begin:ll_end) = phi(ij,ll_begin:ll_end)-phis(ij)
568          ulon_phy(offset,ll_begin:ll_end) = ulon(ij,ll_begin:ll_end)
569          ulat_phy(offset,ll_begin:ll_end) = ulat(ij,ll_begin:ll_end)
570          q_phy(offset,ll_begin:ll_end,:) = q(ij,ll_begin:ll_end,:)
571          wflux_phy(offset,ll_begin:ll_end) = wflux(ij,ll_begin:ll_end)
572        ENDDO
573      ENDDO
574     
575     END SUBROUTINE grid_icosa_to_physics
576
577  END SUBROUTINE physics
578 
579 
580END MODULE physics_lmdz_generic_mod
581   
582 
Note: See TracBrowser for help on using the repository browser.