Ignore:
Timestamp:
06/22/16 15:43:40 (8 years ago)
Author:
dumas
Message:

OpenMP parallelization in Temperature-routines, deformation_mod_2lois and resol_adv_diff_2D-sept2009

File:
1 edited

Legend:

Unmodified
Added
Removed
  • trunk/SOURCES/Temperature-routines/temp_col.f90

    r24 r73  
    1414!> 
    1515 
    16 Subroutine  Temp_col(Ii,Jj,Newtempcol_m,Ct_m,Flot_m,H_m,Tbmer_m) 
    17  
    18   Use Icetemp_declar 
     16Subroutine  Temp_col(Nz,Nzm,Nn,Newtempcol,Ctij,Flotij,Hij,Tbmerij,Ibase,Tij,Tpmpij,Aa,Bb,Cc,Rr,Hh,Ncond,Dee,Cm, & 
     17                                Dzm,Phidij,Ghfij,Ifail1,Tbdotij,DTT,Dou,Tsij,Bmeltij) 
     18 
     19!  Use Icetemp_declar 
    1920  Use Tridiagmod 
    2021   
    2122   
    2223  Implicit None 
    23   integer :: Ii,Jj 
    24   Real,Dimension(Nx,Ny,Nz+Nzm),Intent(Out)  :: Newtempcol_m        !< Tableau 1d Vert. Des Temperatures En Sortie 
    25   Real,Dimension(Nz),Intent(In)             :: Ct_m                !< Tableau 1d Vert. Conductivite Thermique 
    26   Logical,Intent(In)                        :: Flot_m              !< Vrai Si Flot_mtant (Test D'Archimede) 'O' 
    27   Real,Intent(In)                           :: H_m                 !< Ice Thickness  'O' 
    28   Real,Intent(In)                           :: Tbmer_m             !< Temperature De La Mer A La Base De L'Ice Shelf 
    29  
    30  
    31   If (H_m.Gt.10.) Then 
     24  integer,intent(in) :: Nz,Nzm,Nn 
     25  Real,Dimension(Nz+Nzm),Intent(Out)  :: Newtempcol        !< Tableau 1d Vert. Des Temperatures En Sortie 
     26  Real,Dimension(Nz),Intent(In)       :: Ctij                !< Tableau 1d Vert. Conductivite Thermique 
     27  Logical,Intent(In)                  :: Flotij              !< Vrai Si Flotijtant (Test D'Archimede) 'O' 
     28  Real,Intent(In)                     :: Hij                 !< Ice Thickness  'O' 
     29  Real,Intent(In)                     :: Tbmerij             !< Temperature De La Mer A La Base De L'Ice Shelf 
     30  integer,intent(inout)               :: Ibase 
     31  real,dimension(Nz+Nzm),intent(inout) :: Tij             !< Temperature sur la colonne 
     32  real,dimension(Nz),intent(in)       :: Tpmpij               !< pressure melting point temperature in ice sheet 'o' 
     33  Real,Dimension(Nn),intent(inout) :: Aa    !< Work Arrays For Tridiag  !Dim Nn 
     34  Real,Dimension(Nn),intent(inout) :: Bb    !< Work Arrays For Tridiag  !Dim Nn 
     35  Real,Dimension(Nn),intent(inout) :: Cc    !< Work Arrays For Tridiag  !Dim Nn 
     36  Real,Dimension(Nn),intent(inout) :: Rr    !< Work Arrays For Tridiag  !Dim Nn 
     37  Real,Dimension(Nn),intent(inout) :: Hh    !< Work Arrays For Tridiag  !Dim Nn 
     38  Integer,intent(in)               :: Ncond !< Switch : Ncond=1 Thermal Conductivity In Mantle 
     39  real, intent(in) :: Dee 
     40  real,intent(in) :: Cm 
     41  real,intent(in) :: Dzm 
     42  real,intent(inout) :: Phidij     !< flux de chaleur lie a la deformation et glissement basal 
     43  real,intent(in) :: Ghfij      !< flux geothermique 
     44  integer,intent(out) :: Ifail1 !< Diagnostique erreur dans Tridiag 
     45  real, intent(out) :: Tbdotij  !< variation in time of basal temperature 
     46  real, intent(in) :: DTT 
     47  Real, intent(inout) :: Dou 
     48  real, intent(in) :: Tsij        !< surface ice temperature  'o' 
     49  real, intent(out) :: Bmeltij    !< fonte basale 
     50   
     51   
     52   
     53  ! variables locales 
     54  integer :: K 
     55  Real :: Dzi,Tss 
     56  Real,Dimension(Nzm+1) :: Abis,Bbis,Cbis,Rbis,Hbis    
     57 
     58  If (Hij.Gt.10.) Then 
    3259 
    3360     ! Conditions Aux Limites a La Base De La Glace 
     
    3562 
    3663     ! Base Froide  
    37      If ( ((Ibase(Ii,Jj).Eq.1).Or.(Ibase(Ii,Jj).Eq.4) & 
    38           .Or. ((Ibase(Ii,Jj).Eq.5).And.(T(Ii,Jj,Nz).Lt.Tpmp(Ii,Jj,Nz)))) & 
    39           .And..Not.Flot_m) Then 
    40  
    41         Ibase(Ii,Jj)=1 
     64     If ( ((Ibase.Eq.1).Or.(Ibase.Eq.4) & 
     65          .Or. ((Ibase.Eq.5).And.(Tij(Nz).Lt.Tpmpij(Nz)))) & 
     66          .And..Not.Flotij) Then 
     67 
     68        Ibase=1 
    4269        Bb(Nz)=1. 
    4370 
    4471        If (Ncond.Eq.1) Then ! Avec Socle 
    45            Dzi=H_m*Dee*Cm/Ct_m(Nz) 
     72           Dzi=Hij*Dee*Cm/Ctij(Nz) 
    4673           Aa(Nz)=-Dzm/(Dzm+Dzi) 
    4774           Cc(Nz)=-Dzi/(Dzm+Dzi) 
    48            Rr(Nz)=H_m*Dee*Phid(Ii,Jj)*Dzm/Ct_m(Nz)/(Dzm+Dzi) 
     75           Rr(Nz)=Hij*Dee*Phidij*Dzm/Ctij(Nz)/(Dzm+Dzi) 
    4976        Else                 ! Sans Socle 
    5077           Aa(Nz)=-1. 
    5178           Cc(Nz)=0. 
    52            Rr(Nz)=-(Ghf(Ii,Jj)-Phid(Ii,Jj))/Ct_m(Nz)*H_m*Dee 
     79           Rr(Nz)=-(Ghfij-Phidij)/Ctij(Nz)*Hij*Dee 
    5380        Endif 
    5481 
     
    5784        ! Base Temperee Ou Shelf 
    5885        ! ----------------------                
    59         If (Ibase(Ii,Jj).Eq.5)  Ibase(Ii,Jj)=2 
    60         Ibase(Ii,Jj)=Max(Ibase(Ii,Jj),2) 
     86        If (Ibase.Eq.5)  Ibase=2 
     87        Ibase=Max(Ibase,2) 
    6188        Aa(Nz)=0. 
    6289        Bb(Nz)=1. 
    6390        Cc(Nz)=0. 
    6491 
    65         If (.Not.Flot_m) Then 
    66            Rr(Nz)=Tpmp(Ii,Jj,Nz) 
    67         Else 
    68            Rr(Nz)=Tbmer_m 
     92        If (.Not.Flotij) Then 
     93           Rr(Nz)=Tpmpij(Nz) 
     94        Else 
     95           Rr(Nz)=Tbmerij 
    6996        End If 
    7097 
     
    74101     If (Ncond.Eq.1) Then ! Avec Socle 
    75102        Do K=Nz+1,nz+nzm-1 
    76            Rr(K)=T(Ii,Jj,K) 
     103           Rr(K)=Tij(K) 
    77104        End Do 
    78105        Call Tridiag (Aa,Bb,Cc,Rr,Hh,nz+nzm,Ifail1) 
     
    87114     If (Ncond.Eq.0) Then 
    88115        Do K=Nz+1,nz+nzm 
    89            Hh(K)=Hh(Nz)-Dzm*(K-Nz)*Ghf(Ii,Jj)/Cm 
     116           Hh(K)=Hh(Nz)-Dzm*(K-Nz)*Ghfij/Cm 
    90117        End Do 
    91118     Endif 
    92119 
    93120     Do K=1,nz+nzm 
    94         Newtempcol_m(Ii,Jj,K)=Hh(K) 
     121        Newtempcol(K)=Hh(K) 
    95122     End Do 
    96      Tbdot(Ii,Jj)=(Newtempcol_m(Ii,Jj,Nz)-T(Ii,Jj,Nz))/Dtt 
     123     Tbdotij=(Newtempcol(Nz)-Tij(Nz))/Dtt 
    97124 
    98125     !      ------------------------------------------------- Glace Tres Fine (H<10m Ou H=0) 
    99   Else If (H_m.Le.10.) Then 
     126  Else If (Hij.Le.10.) Then 
    100127 
    101128     !       Pour Eviter Des Problemes Lorsque H Est Inferieur A 10 M 
     
    104131 
    105132     ! ........................................ Avec Calcul Dans Le Socle 
    106      If (Ncond.Eq.1.And..Not.Flot_m) Then    
    107         If (H_m.Gt.0.) Then 
     133     If (Ncond.Eq.1.And..Not.Flotij) Then    
     134        If (Hij.Gt.0.) Then 
    108135           ! Gradient Dans Le Socle 
    109            Dou=(T(Ii,Jj,Nz+1)-T(Ii,Jj,Nz))/Dzm*Cm      
    110            Dou=Dou/Ct_m(Nz)*Dee*H_m 
    111  
    112            Tss=Min(0.,Ts(Ii,Jj)) 
    113            Do K=1,Nz 
    114               Newtempcol_m(Ii,Jj,K)=Tss+Dou*(K-1.) 
     136           Dou=(Tij(Nz+1)-Tij(Nz))/Dzm*Cm      
     137           Dou=Dou/Ctij(Nz)*Dee*Hij 
     138 
     139           Tss=Min(0.,Tsij) 
     140           Do K=1,Nz 
     141              Newtempcol(K)=Tss+Dou*(K-1.) 
    115142           End Do 
    116143 
    117144        Else 
    118145           ! Pas De Glace 
    119            Tss=Ts(Ii,Jj) 
    120            Do K=1,Nz 
    121               Newtempcol_m(Ii,Jj,K)=Tss 
     146           Tss=Tsij 
     147           Do K=1,Nz 
     148              Newtempcol(K)=Tss 
    122149           End Do 
    123150        End If 
    124151 
    125152        ! Calcul Dans Le Socle Meme S'Il N'Y A Pas De Glace 
    126         Rr(Nz)=Newtempcol_m(Ii,Jj,Nz) 
     153        Rr(Nz)=Newtempcol(Nz) 
    127154        Aa(Nz)=0 
    128155        Bb(Nz)=1 
     
    130157 
    131158        Do K=Nz+1,nz+nzm-1 
    132            Rr(K)=T(Ii,Jj,K) 
     159           Rr(K)=Tij(K) 
    133160        End Do 
    134161 
     
    150177     Else 
    151178        ! Calotte Posee                     
    152         If ((H_m.Gt.0.).And..Not.Flot_m) Then 
    153            Dou=-Ghf(Ii,Jj)/Ct_m(Nz)*Dee*H_m 
    154            Tss=Min(0.,Ts(Ii,Jj)) 
    155            Do K=1,Nz 
    156               Newtempcol_m(Ii,Jj,K)=Tss+Dou*(K-1.) 
     179        If ((Hij.Gt.0.).And..Not.Flotij) Then 
     180           Dou=-Ghfij/Ctij(Nz)*Dee*Hij 
     181           Tss=Min(0.,Tsij) 
     182           Do K=1,Nz 
     183              Newtempcol(K)=Tss+Dou*(K-1.) 
    157184           End Do 
    158185 
    159186           !  Shelf 
    160         Else If ((H_m.Gt.0.).And.Flot_m) Then 
    161            Tss=Min(0.,Ts(Ii,Jj)) 
    162            Dou=(Tbmer_m-Tss)*Dee 
    163            Do K=1,Nz 
    164               Newtempcol_m(Ii,Jj,K)=Tss+Dou*(K-1.) 
     187        Else If ((Hij.Gt.0.).And.Flotij) Then 
     188           Tss=Min(0.,Tsij) 
     189           Dou=(Tbmerij-Tss)*Dee 
     190           Do K=1,Nz 
     191              Newtempcol(K)=Tss+Dou*(K-1.) 
    165192           End Do 
    166193 
    167194           ! Pas De Glace 
    168195        Else 
    169            Tss=Ts(Ii,Jj)       
    170            Do K=1,Nz 
    171               Newtempcol_m(Ii,Jj,K)=Tss 
     196           Tss=Tsij       
     197           Do K=1,Nz 
     198              Newtempcol(K)=Tss 
    172199           End Do 
    173200        End If 
    174201 
    175202        ! Temperature Dans Le Socle, Lineaire Avec Le Gradient Ghf 
    176         If (.Not.Flot_m) Then 
     203        If (.Not.Flotij) Then 
    177204           Do K=Nz+1,nz+nzm 
    178               Hh(K)=Newtempcol_m(Ii,Jj,Nz)-Dzm*(K-Nz)*Ghf(Ii,Jj)/Cm 
     205              Hh(K)=Newtempcol(Nz)-Dzm*(K-Nz)*Ghfij/Cm 
    179206           End Do 
    180207        Else 
    181208           Do K=Nz+1,nz+nzm 
    182               Hh(K)=Tbmer_m-Dzm*(K-Nz)*Ghf(Ii,Jj)/Cm 
     209              Hh(K)=Tbmerij-Dzm*(K-Nz)*Ghfij/Cm 
    183210           End Do 
    184211        End If 
     
    189216 
    190217     Do K=Nz+1,nz+nzm 
    191         Newtempcol_m(Ii,Jj,K)=Hh(K) 
     218        Newtempcol(K)=Hh(K) 
    192219     End Do 
    193220 
    194      Tbdot(Ii,Jj)=(Newtempcol_m(Ii,Jj,Nz)-T(Ii,Jj,Nz))/Dtt 
    195      Bmelt(Ii,Jj)=0. 
    196      Ibase(Ii,Jj)=5 
    197      Phid(Ii,Jj)=0. 
     221     Tbdotij=(Newtempcol(Nz)-Tij(Nz))/Dtt 
     222     Bmeltij=0. 
     223     Ibase=5 
     224     Phidij=0. 
    198225 
    199226 
Note: See TracChangeset for help on using the changeset viewer.