source: trunk/LMDZ.COMMON/libf/evolution/atmosphere.F90 @ 4076

Last change on this file since 4076 was 4074, checked in by jbclement, 11 days ago

PEM:

  • Correct management of H2O ice tendency in 1D when there is not enough ice anymore.
  • Clean initialization of allocatable module arrays (especially needed when no slope)
  • One more renaming for consistency + few small updates thoughout the code.

JBC

File size: 20.0 KB
Line 
1MODULE atmosphere
2!-----------------------------------------------------------------------
3! NAME
4!     atmosphere
5!
6! DESCRIPTION
7!     Contains global parameters used for the atmosphere.
8!
9! AUTHORS & DATE
10!     JB Clement, 12/2025
11!
12! NOTES
13!
14!-----------------------------------------------------------------------
15
16! DEPENDENCIES
17! ------------
18use numerics, only: dp, di, k4
19
20! DECLARATION
21! -----------
22implicit none
23
24! PARAMETERS
25! ----------
26character(9),                          parameter :: z2sigdef_name = 'z2sig.def'
27logical(k4),                           protected :: hybrid_alt_coord   ! Flag fo hybrid coordinates
28real(dp), dimension(:),   allocatable, protected :: ap             ! Hybrid (pressure contribution) coordinate at layer interfaces [Pa]
29real(dp), dimension(:),   allocatable, protected :: bp             ! Hybrid (sigma contribution) coordinate at layer interfaces
30real(dp), dimension(:),   allocatable, protected :: ps_PCM         ! Surface pressure in the "start.nc" at the beginning [Pa]
31real(dp), dimension(:,:), allocatable, protected :: teta_PCM       ! Potential temperature in the "start.nc" at the beginning [K]
32real(dp), dimension(:,:), allocatable, protected :: u_PCM          ! Zonal wind [m/s]
33real(dp), dimension(:,:), allocatable, protected :: v_PCM          ! Meridional wind [m/s]
34real(dp)                                         :: CO2cond_ps_PCM ! Coefficient to control the surface pressure change
35
36contains
37!+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
38
39!=======================================================================
40SUBROUTINE ini_atmosphere()
41!-----------------------------------------------------------------------
42! NAME
43!     ini_atmosphere
44!
45! DESCRIPTION
46!     Initialize the parameters of module 'atmosphere'.
47!
48! AUTHORS & DATE
49!     JB Clement, 12/2025
50!
51! NOTES
52!
53!-----------------------------------------------------------------------
54
55! DEPENDENCIES
56! ------------
57use geometry, only: nlayer, ngrid
58
59! DECLARATION
60! -----------
61implicit none
62
63! CODE
64! ----
65if (.not. allocated(ap)) allocate(ap(nlayer + 1))
66if (.not. allocated(bp)) allocate(bp(nlayer + 1))
67if (.not. allocated(ps_PCM)) allocate(ps_PCM(ngrid))
68if (.not. allocated(teta_PCM)) allocate(teta_PCM(ngrid,nlayer))
69if (.not. allocated(u_PCM)) allocate(u_PCM(ngrid,nlayer))
70if (.not. allocated(v_PCM)) allocate(v_PCM(ngrid,nlayer))
71ap(:) = 0._dp
72bp(:) = 0._dp
73ps_PCM(:) = 0._dp
74teta_PCM(:,:) = 0._dp
75u_PCM(:,:) = 0._dp
76v_PCM(:,:) = 0._dp
77
78END SUBROUTINE ini_atmosphere
79!=======================================================================
80
81!=======================================================================
82SUBROUTINE end_atmosphere()
83!-----------------------------------------------------------------------
84! NAME
85!     end_atmosphere
86!
87! DESCRIPTION
88!     Deallocate atmosphere arrays.
89!
90! AUTHORS & DATE
91!     JB Clement, 12/2025
92!
93! NOTES
94!
95!-----------------------------------------------------------------------
96
97! DECLARATION
98! -----------
99implicit none
100
101! CODE
102! ----
103if (allocated(ap)) deallocate(ap)
104if (allocated(bp)) deallocate(bp)
105if (allocated(ps_PCM)) deallocate(ps_PCM)
106if (allocated(teta_PCM)) deallocate(teta_PCM)
107if (allocated(u_PCM)) deallocate(u_PCM)
108if (allocated(v_PCM)) deallocate(v_PCM)
109
110END SUBROUTINE end_atmosphere
111!=======================================================================
112
113!=======================================================================
114SUBROUTINE set_atmosphere_config(hybrid_alt_coord_in)
115!-----------------------------------------------------------------------
116! NAME
117!     set_atmosphere_config
118!
119! DESCRIPTION
120!     Setter for 'atmosphere' configuration parameters.
121!
122! AUTHORS & DATE
123!     JB Clement, 02/2026
124!
125! NOTES
126!
127!-----------------------------------------------------------------------
128
129! DEPENDENCIES
130! ------------
131use utility, only: bool2str
132use display, only: print_msg
133
134! DECLARATION
135! -----------
136implicit none
137
138! ARGUMENTS
139! ---------
140logical(k4), intent(in) :: hybrid_alt_coord_in
141
142! CODE
143! ----
144hybrid_alt_coord = hybrid_alt_coord_in
145call print_msg('hybrid_alt_coord = '//bool2str(hybrid_alt_coord))
146
147END SUBROUTINE set_atmosphere_config
148!=======================================================================
149
150!=======================================================================
151SUBROUTINE set_ap(ap_in)
152!-----------------------------------------------------------------------
153! NAME
154!     set_ap
155!
156! DESCRIPTION
157!     Setter for 'ap'.
158!
159! AUTHORS & DATE
160!     JB Clement, 12/2025
161!
162! NOTES
163!
164!-----------------------------------------------------------------------
165
166! DECLARATION
167! -----------
168implicit none
169
170! ARGUMENTS
171! ---------
172real(dp), dimension(:), intent(in) :: ap_in
173
174! CODE
175! ----
176ap(:) = ap_in(:)
177
178END SUBROUTINE set_ap
179!=======================================================================
180
181!=======================================================================
182SUBROUTINE set_bp(bp_in)
183!-----------------------------------------------------------------------
184! NAME
185!     set_bp
186!
187! DESCRIPTION
188!     Setter for 'bp'.
189!
190! AUTHORS & DATE
191!     JB Clement, 12/2025
192!
193! NOTES
194!
195!-----------------------------------------------------------------------
196
197! DECLARATION
198! -----------
199implicit none
200
201! ARGUMENTS
202! ---------
203real(dp), dimension(:), intent(in) :: bp_in
204
205! CODE
206! ----
207bp(:) = bp_in(:)
208
209END SUBROUTINE set_bp
210!=======================================================================
211
212!=======================================================================
213SUBROUTINE set_ps_PCM(ps_PCM_in)
214!-----------------------------------------------------------------------
215! NAME
216!     set_ps_PCM
217!
218! DESCRIPTION
219!     Setter for 'ps_PCM'.
220!
221! AUTHORS & DATE
222!     JB Clement, 12/2025
223!
224! NOTES
225!
226!-----------------------------------------------------------------------
227
228! DEPENDENCIES
229! ------------
230use geometry, only: dyngrd2vect, nlon, nlat, ngrid
231
232! DECLARATION
233! -----------
234implicit none
235
236! ARGUMENTS
237! ---------
238real(dp), dimension(:,:), intent(in) :: ps_PCM_in
239
240! CODE
241! ----
242call dyngrd2vect(nlon + 1,nlat,ngrid,ps_PCM_in,ps_PCM)
243
244END SUBROUTINE set_ps_PCM
245!=======================================================================
246
247!=======================================================================
248SUBROUTINE set_teta_PCM(teta_PCM_in)
249!-----------------------------------------------------------------------
250! NAME
251!     set_teta_PCM
252!
253! DESCRIPTION
254!     Setter for 'teta_PCM'.
255!
256! AUTHORS & DATE
257!     JB Clement, 01/2026
258!
259! NOTES
260!
261!-----------------------------------------------------------------------
262
263! DEPENDENCIES
264! ------------
265use geometry, only: dyngrd2vect, nlon, nlat, ngrid, nlayer
266
267! DECLARATION
268! -----------
269implicit none
270
271! ARGUMENTS
272! ---------
273real(dp), dimension(:,:,:), intent(in) :: teta_PCM_in
274
275! LOCAL VARIABLES
276! ---------------
277integer(di) :: l
278
279! CODE
280! ----
281do l = 1,nlayer
282    call dyngrd2vect(nlon + 1,nlat,ngrid,teta_PCM_in(:,:,l),teta_PCM(:,l))
283end do
284
285END SUBROUTINE set_teta_PCM
286!=======================================================================
287
288!=======================================================================
289SUBROUTINE set_u_PCM(u_PCM_in)
290!-----------------------------------------------------------------------
291! NAME
292!     set_u_PCM
293!
294! DESCRIPTION
295!     Setter for 'u_PCM'.
296!
297! AUTHORS & DATE
298!     JB Clement, 01/2026
299!
300! NOTES
301!
302!-----------------------------------------------------------------------
303
304! DEPENDENCIES
305! ------------
306use geometry, only: dyngrd2vect, nlon, nlat, ngrid, nlayer
307
308! DECLARATION
309! -----------
310implicit none
311
312! ARGUMENTS
313! ---------
314real(dp), dimension(:,:,:), intent(in) :: u_PCM_in
315
316! LOCAL VARIABLES
317! ---------------
318integer(di) :: l
319
320! CODE
321! ----
322do l = 1,nlayer
323    call dyngrd2vect(nlon + 1,nlat,ngrid,u_PCM_in(:,:,l),u_PCM(:,l))
324end do
325
326END SUBROUTINE set_u_PCM
327!=======================================================================
328
329!=======================================================================
330SUBROUTINE set_v_PCM(v_PCM_in)
331!-----------------------------------------------------------------------
332! NAME
333!     set_v_PCM
334!
335! DESCRIPTION
336!     Setter for 'v_PCM'.
337!
338! AUTHORS & DATE
339!     JB Clement, 01/2026
340!
341! NOTES
342!
343!-----------------------------------------------------------------------
344
345! DEPENDENCIES
346! ------------
347use geometry, only: dyngrd2vect, nlon, nlat, ngrid, nlayer
348
349! DECLARATION
350! -----------
351implicit none
352
353! ARGUMENTS
354! ---------
355real(dp), dimension(:,:,:), intent(in) :: v_PCM_in
356
357! LOCAL VARIABLES
358! ---------------
359integer(di) :: l
360
361! CODE
362! ----
363do l = 1,nlayer
364    call dyngrd2vect(nlon + 1,nlat,ngrid,v_PCM_in(:,:,l),v_PCM(:,l))
365end do
366
367END SUBROUTINE set_v_PCM
368!=======================================================================
369
370!=======================================================================
371SUBROUTINE compute_alt_coord(pa,preff,ap,bp,aps,bps)
372!-----------------------------------------------------------------------
373! NAME
374!     compute_alt_coord
375!
376! DESCRIPTION
377!     Compute values of the coefficients to define the altitude
378!     coordinates.
379!
380! AUTHORS & DATE
381!     JB Clement, 02/2026
382!
383! NOTES
384!     Taken from 'disvert_noterre' written by F. Forget, Y. Wanherdrick
385!     and P. Levan in the dynamics.
386!
387!     'pa' and 'preff' are read from "start.nc"/"start1D.txt".
388!
389!     WARNING: to compute the values at mid_layer, there is an arbitrary
390!     decision here to put the middle of the layer half-way (pressure-wise)
391!     between boundaries.
392!     In addition, for the last layer (the top of which is at zero pressure)
393!     P(llm), we choose to put it with the same interval (in log(pressure))
394!     from P(llm-1), than there is between between P(llm-1) and P(llm-2)
395!     This choice is quite specific and must be compliant with the PCM!
396!-----------------------------------------------------------------------
397
398! DEPENDENCIES
399! ------------
400use geometry, only: nlayer
401use display,  only: print_msg
402use stoppage, only: stop_clean
403use utility,  only: real2str
404
405! DECLARATION
406! -----------
407implicit none
408
409! ARGUMENTS
410! ---------
411real(dp)                                               :: pa, preff
412real(dp), dimension(nlayer + 1), intent(out)           :: ap, bp
413real(dp), dimension(nlayer),     intent(out), optional :: aps, bps
414
415! LOCAL VAIRABLES
416! ---------------
417integer(di)                     :: i
418real(dp)                        :: scaleheight, newsig
419real(dp), dimension(nlayer + 1) :: sig
420
421! CODE
422! ----
423! Read vertical grid definition
424call read_z2sig(scaleheight,sig)
425call print_msg('Pa    [Pa] = '//real2str(pa))
426call print_msg('Preff [Pa] = '//real2str(preff))
427
428! Compute ap and bp coefficients
429if (hybrid_alt_coord) then
430    call print_msg('> Defining hybrid altitude coordinates')
431    do i = 1,nlayer
432        call compute_hybrid_sig(sig(i),pa,preff,newsig)
433        bp(i) = exp(1._dp - 1._dp/(newsig**2))
434        ap(i) = pa*(newsig - bp(i))
435    end do
436    ap(nlayer + 1) = 0._dp
437    bp(nlayer + 1) = 0._dp
438else
439    call print_msg('> Defining sigma altitude coordinates')
440    ap(:) = 0._dp
441    bp(1:nlayer) = sig(1:nlayer)
442    bp(nlayer + 1) = 0._dp
443end if
444
445! Compute aps and bps coefficients (mid-layer) if asked
446if (present(aps) .neqv. present(bps)) call stop_clean(__FILE__,__LINE__,"'aps' and 'bps' must be computed together!",1)
447if (.not. present(aps)) return
448
449do i = 1,nlayer - 1
450    aps(i) =  0.5_dp*(ap(i) + ap(i + 1))
451    bps(i) =  0.5_dp*(bp(i) + bp(i + 1))
452end do
453if (hybrid_alt_coord) then
454    aps(nlayer) = aps(nlayer - 1)**2/aps(nlayer - 2)
455    bps(nlayer) = 0.5_dp*(bp(nlayer) + bp(nlayer + 1))
456else
457    aps(nlayer) = 0._dp
458    bps(nlayer) = bps(nlayer - 1)**2/bps(nlayer - 2)
459end if
460
461END SUBROUTINE compute_alt_coord
462!=======================================================================
463
464!=======================================================================
465SUBROUTINE compute_hybrid_sig(sig,pa,preff,newsig)
466!-----------------------------------------------------------------------
467! NAME
468!     compute_hybrid_sig
469!
470! DESCRIPTION
471!     Compute modified sigma value to keep vertical coordinates described
472!     in "z2sig.def" when defining hybrid coordinates.
473!
474! AUTHORS & DATE
475!     JB Clement, 01/2026
476!
477! NOTES
478!     Taken from 'sig_hybrid' written by F. Forget (2002) in the
479!     dynamics.
480!
481!     Knowing sig (the "sigma" levels where we want to put the layers),
482!     it solves 'newsig' such that:
483!       (1 - pa/preff)*exp(1 - 1/newsig^2) + (pa/preff)*newsig = sig
484!     which cannot be done analyticaly.
485!     => needs for an iterative procedure
486!
487!     Information: when exp(1 - 1/x**2) becomes << x
488!           x      exp(1 - 1/x**2)/x
489!           1           1
490!           0.68       0.5
491!           0.5        1.E-1
492!           0.391      1.E-2
493!           0.333      1.E-3
494!           0.295      1.E-4
495!           0.269      1.E-5
496!           0.248      1.E-6
497!     => one can thus use newsig = sig*preff/pa if sig*preff/pa < 0.25
498!-----------------------------------------------------------------------
499
500! DEPENDENCIES
501! ------------
502use stoppage, only: stop_clean
503
504! DECLARATION
505! -----------
506implicit none
507
508! ARGUMENTS
509! ---------
510real(dp), intent(in)  :: sig, pa, preff
511real(dp), intent(out) :: newsig
512
513! LOCAL VARIABLES
514! ---------------
515real(dp)    :: x1, x2, F
516integer(di) :: j
517
518! CODE
519! ----
520newsig = sig
521x1 = 0._dp
522x2 = 1._dp
523
524if (sig >= 1.) then
525    newsig = sig
526else if (sig*preff/pa >= 0.25_dp) then
527    do j = 1,9999 ! Overkill maximum number of iterations
528        F = ((1._dp - pa/preff)*exp(1._dp - 1._dp/newsig**2) + (pa/preff)*newsig)/sig
529
530        if (F > 1.) then
531            x2 = newsig
532            newsig = 0.5_dp*(x1 + newsig)
533        else
534        x1 = newsig
535        newsig = 0.5_dp*(x2 + newsig)
536        end if
537        ! Convergence is assumed if sig is within 0.01m (in pseudo-altitude) of the target value
538        if (abs(10._dp * log(F)) < 1.e-5_dp) exit
539    end do
540else !if (sig*preff/pa <= 0.25_dp) then
541    newsig = sig*preff/pa
542end if
543
544END SUBROUTINE compute_hybrid_sig
545!=======================================================================
546
547!=======================================================================
548SUBROUTINE read_z2sig(scaleheight,sig)
549!-----------------------------------------------------------------------
550! NAME
551!     read_z2sig
552!
553! DESCRIPTION
554!     Read vertical grid definition from z2sig.def and compute sigma
555!     levels.
556!
557! AUTHORS & DATE
558!     JB Clement, 01/2026
559!
560! NOTES
561!     Adapted from 'disvert_noterre' written by F. Forget, Y.
562!     Wanherdrick and P. Levan in the dynamics.
563!-----------------------------------------------------------------------
564
565! DEPENDENCIES
566! ------------
567use geometry, only: nlayer
568use stoppage, only: stop_clean
569use display,  only: print_msg
570use utility,  only: real2str
571
572! DECLARATION
573! -----------
574implicit none
575
576! ARGUMENTS
577! ---------
578real(dp),                        intent(out) :: scaleheight
579real(dp), dimension(nlayer + 1), intent(out) :: sig
580
581! LOCAL VARIABLES
582! ---------------
583real(dp), dimension(nlayer) :: zsig
584integer(di)                 :: l, ierr, funit
585logical(k4)                 :: here
586
587! CODE
588! ----
589! Open "z2sig.def"
590inquire(file = z2sigdef_name,exist = here)
591if (.not. here) call stop_clean(__FILE__,__LINE__,'cannot find required file "'//z2sigdef_name//'"!',1)
592call print_msg('> Reading "'//z2sigdef_name//'"')
593open(newunit = funit,file = z2sigdef_name,status = "old",action = "read",iostat = ierr)
594if (ierr /= 0) call stop_clean(__FILE__,__LINE__,'error opening "'//z2sigdef_name//'"',ierr)
595
596! Read scale height
597read(funit,*,iostat = ierr) scaleheight
598if (ierr /= 0) call stop_clean(__FILE__,__LINE__,"error reading 'scaleheight'!",1)
599
600! Read pseudo-altitudes
601do l = 1,nlayer
602    read(funit,*,iostat = ierr) zsig(l)
603    if (ierr /= 0) call stop_clean(__FILE__,__LINE__,'"'//z2sigdef_name//'" too short for the number of altitude levels!',1)
604end do
605
606! Check for extra lines in the file (warning only)
607read(funit,*,iostat = ierr)
608if (ierr == 0) call print_msg('Warning: "'//z2sigdef_name//'" has more lines than expected!')
609
610! Compute sigma values
611sig(1) = 1.
612do l = 2,nlayer
613    sig(l) = 0.5_dp*(exp(-zsig(l)/scaleheight) + exp(-zsig(l - 1)/scaleheight))
614end do
615sig(nlayer + 1) = 0._dp
616
617! Close
618close(funit)
619call print_msg('Scale height [km] = '//real2str(scaleheight))
620
621END SUBROUTINE read_z2sig
622!=======================================================================
623
624!=======================================================================
625SUBROUTINE build4PCM_atmosphere(ps_avg,ps_dev,ps_avg_glob,ps_avg_glob_ini,ps4PCM,pa4PCM,preff4PCM,teta4PCM,air_mass4PCM)
626!-----------------------------------------------------------------------
627! NAME
628!     build4PCM_atmosphere
629!
630! DESCRIPTION
631!     Reconstructs the pressure, potential temperature and tha air mass
632!     for the PCM.
633!
634! AUTHORS & DATE
635!     JB Clement, 01/2026
636!
637! NOTES
638!     Air mass computation taken from 'massdair' and 'pression' written
639!     by P. Le Van and F. Hourdin in the dynamics.
640!-----------------------------------------------------------------------
641
642! DEPENDENCIES
643! ------------
644use geometry, only: ngrid, nlayer, cell_area, vect2dyngrd
645use physics,  only: g, rcp
646use display,  only: print_msg
647
648! DECLARATION
649! -----------
650implicit none
651
652! ARGUMENTS
653! ---------
654real(dp), dimension(:),   intent(in)  :: ps_avg, ps_dev
655real(dp),                 intent(in)  :: ps_avg_glob, ps_avg_glob_ini
656real(dp), dimension(:,:), intent(out) :: air_mass4PCM, teta4PCM
657real(dp), dimension(:),   intent(out) :: ps4PCM
658real(dp),                 intent(out) :: pa4PCM, preff4PCM
659
660! LOCAL VARIABLES
661! ---------------
662real(dp), dimension(ngrid,nlayer + 1) :: p ! 3D pressure field
663integer(di)                           :: l
664
665! CODE
666! ----
667call print_msg('> Building pressure for the PCM')
668! The pressure deviation is rescaled to avoid disproportionate oscillations in case of huge change of average pressure during the PEM run
669ps4PCM(:) = ps_avg(:) + ps_dev(:)*ps_avg_glob/ps_avg_glob_ini
670pa4PCM = ps_avg_glob_ini/30. ! For now the altitude grid is not changed
671preff4PCM = ps_avg_glob_ini  ! For now the altitude grid is not changed
672
673! Correction on teta due to surface pressure changes
674call print_msg('> Building potential temperature for the PCM')
675do l = 1,nlayer
676    teta4PCM(:,l) = teta_PCM(:,l)*ps4PCM(:)**rcp
677end do
678
679! Compute atmospheric pressure
680call print_msg('> Building air mass for the PCM')
681do l = 1,nlayer + 1
682    p(:,l) = ap(l) + bp(l)*ps4PCM(:)
683end do
684
685! Compute air mass
686do l = 1,nlayer
687    air_mass4PCM(:,l) = cell_area(:)*(p(:,l) - p(:,l + 1))/g
688end do
689
690END SUBROUTINE build4PCM_atmosphere
691!=======================================================================
692
693!=======================================================================
694SUBROUTINE evolve_pressure(d_co2ice,delta_co2_ads,do_sorption,ps_avg_glob_old,ps_avg_glob,ps_avg)
695!-----------------------------------------------------------------------
696! NAME
697!     evolve_pressure
698!
699! DESCRIPTION
700!     Evolve pressure according to tendency.
701!
702! AUTHORS & DATE
703!     JB Clement, 01/2026
704!
705! NOTES
706!
707!-----------------------------------------------------------------------
708
709! DEPENDENCIES
710! ------------
711use geometry, only: ngrid, nslope, cell_area, total_surface
712use slopes,   only: subslope_dist, def_slope_mean
713use physics,  only: g
714use maths,    only: pi
715use display,  only: print_msg
716use utility,  only: real2str
717
718! DECLARATION
719! -----------
720implicit none
721
722! ARGUMENTS
723! ---------
724real(dp), dimension(:,:), intent(in)    :: d_co2ice
725real(dp), dimension(:),   intent(in)    :: delta_co2_ads
726logical(k4),              intent(in)    :: do_sorption
727real(dp),                 intent(inout) :: ps_avg_glob_old, ps_avg_glob
728real(dp), dimension(:),   intent(inout) :: ps_avg
729
730! LOCAL VARIABLES
731! ---------------
732integer(di) :: i, islope
733
734! CODE
735! ----
736call print_msg("> Evolving the surface pressure")
737ps_avg_glob_old = ps_avg_glob
738do i = 1,ngrid
739    do islope = 1,nslope
740        ps_avg_glob = ps_avg_glob - CO2cond_ps_PCM*g*cell_area(i)*d_co2ice(i,islope)*subslope_dist(i,islope)/cos(pi*def_slope_mean(islope)/180.)/total_surface
741    end do
742    if (do_sorption) ps_avg_glob = ps_avg_glob - g*cell_area(i)*delta_co2_ads(i)/total_surface
743end do
744call print_msg('Global average pressure (old|new) [Pa] = '//real2str(ps_avg_glob_old)//' | '//real2str(ps_avg_glob))
745ps_avg(:) = ps_avg(:)*ps_avg_glob/ps_avg_glob_old
746
747END SUBROUTINE evolve_pressure
748!=======================================================================
749
750END MODULE atmosphere
Note: See TracBrowser for help on using the repository browser.