[4] | 1 | !> \file conserv-mass-adv-diff_sept2009_mod.f90 |
---|
| 2 | !! Module qui calcule l evolution de l epaisseur en resolvant |
---|
| 3 | !! une equation d'advection diffusion |
---|
| 4 | !< |
---|
| 5 | |
---|
| 6 | !> \namespace equat_adv_diff_2D_vect |
---|
| 7 | !! Calcule l evolution de l epaisseur en resolvant |
---|
| 8 | !! une equation d'advection diffusion |
---|
| 9 | !! @note Version avec : call resol_adv_diff_2D_vect (Dmx,Dmy,advmx,advmy,H_p,i_Hp,bilmass,vieuxH,H) |
---|
| 10 | !! @note - le terme vecteur est passe en dummy parametres |
---|
| 11 | !! l epaisseur peut etre imposee |
---|
| 12 | !! \author CatRitz |
---|
| 13 | !! \date june 2009 |
---|
| 14 | !! @note Used modules: |
---|
| 15 | !! @note - use module3D_phy |
---|
| 16 | !! @note - use reso_adv_diff_2D |
---|
| 17 | !! @todo essayer dans le futur la methode de Richard dUX/dx = H dU/dx + U dH/dx |
---|
| 18 | !< |
---|
| 19 | module equat_adv_diff_2D_vect ! Cat nouvelle mouture juin 2009 |
---|
| 20 | use module3D_phy |
---|
| 21 | use reso_adv_diff_2D_vect |
---|
| 22 | use prescribe_H |
---|
| 23 | |
---|
| 24 | implicit none |
---|
| 25 | |
---|
| 26 | real, dimension(nx,ny) :: advmx !< partie advective en x |
---|
| 27 | real, dimension(nx,ny) :: advmy !< partie advective en y |
---|
| 28 | real, dimension(nx,ny) :: dmx !< partie advective en x |
---|
| 29 | real, dimension(nx,ny) :: dmy !< partie advective en y |
---|
| 30 | real, dimension(nx,ny) :: VieuxH !< ancienne valeur de H |
---|
| 31 | real, dimension(nx,ny) :: bilmass !< bilan de masse pour la colonne |
---|
| 32 | logical, dimension(nx,ny) :: zonemx !< pour separer advection-diffusion |
---|
| 33 | logical, dimension(nx,ny) :: zonemy !< pour separer advection-diffusion |
---|
| 34 | real :: adv_frac !< fraction du flux traitee en advection |
---|
| 35 | integer :: itesti |
---|
| 36 | integer :: itour |
---|
| 37 | |
---|
| 38 | |
---|
| 39 | contains |
---|
| 40 | |
---|
| 41 | !> SUBROUTINE: init_icethick |
---|
| 42 | !! Definition et initialisation des parametres qui gerent la repartition advection diffusion |
---|
| 43 | !! @return adv_frac |
---|
| 44 | |
---|
| 45 | subroutine init_icethick |
---|
| 46 | |
---|
| 47 | implicit none |
---|
| 48 | |
---|
| 49 | namelist/mass_conserv/adv_frac,V_limit |
---|
| 50 | |
---|
| 51 | |
---|
| 52 | rewind(num_param) ! pour revenir au debut du fichier param_list.dat |
---|
| 53 | read(num_param,mass_conserv) |
---|
| 54 | |
---|
| 55 | |
---|
| 56 | write(num_rep_42,*)'! Conservation de la masse avec equation advection-diffusion ' |
---|
| 57 | write(num_rep_42,mass_conserv) |
---|
| 58 | write(num_rep_42,*)'! la repartition depend de adv_frac' |
---|
| 59 | write(num_rep_42,*)'! >1 -> advection seule' |
---|
| 60 | write(num_rep_42,*)'! 0 -> diffusion seule' |
---|
| 61 | write(num_rep_42,*)'! 0<*<1 -> fraction de l advection' |
---|
| 62 | write(num_rep_42,*)'! -1 -> zones diffusion + zones advection' |
---|
| 63 | write(num_rep_42,*)'! V_limit depend de la calotte : ' |
---|
| 64 | write(num_rep_42,*)'! typiquement 3000 m/an en Antarctique, 10000 m/an au Groenland' |
---|
| 65 | write(num_rep_42,*)'!___________________________________________________________' |
---|
| 66 | |
---|
| 67 | |
---|
| 68 | call init_reso_adv_diff_2D |
---|
| 69 | call init_prescribe_H ! initialize grounding line mask |
---|
| 70 | |
---|
| 71 | Dx1=1./Dx |
---|
| 72 | itour=0 |
---|
| 73 | |
---|
| 74 | end subroutine init_icethick |
---|
| 75 | |
---|
| 76 | !------------------------------------------------------------------ |
---|
| 77 | !> SUBROUTINE: icethick3 |
---|
| 78 | !! Routine pour le calcule de l'epaisseur de glace |
---|
| 79 | !< |
---|
| 80 | subroutine icethick3 |
---|
| 81 | |
---|
[77] | 82 | !$ USE OMP_LIB |
---|
| 83 | |
---|
[4] | 84 | implicit none |
---|
| 85 | real,dimension(nx,ny) :: Dminx,Dminy |
---|
| 86 | real,dimension(nx,ny) :: Uxdiff,Uydiff ! vitesse due a la diffusion |
---|
| 87 | real aux ! pour le calcul du critere |
---|
| 88 | real maxdia ! sur le pas de temps |
---|
| 89 | integer imaxdia,jmaxdia |
---|
| 90 | |
---|
| 91 | |
---|
| 92 | if (itracebug.eq.1) call tracebug(" Entree dans routine icethick") |
---|
| 93 | |
---|
| 94 | ! Copie de H sur vieuxH |
---|
| 95 | ! -------------------- |
---|
[77] | 96 | !$OMP PARALLEL |
---|
| 97 | !$OMP WORKSHARE |
---|
[4] | 98 | where (mk0(:,:).eq.0) |
---|
| 99 | vieuxH(:,:)=0. |
---|
| 100 | elsewhere (i_Hp(:,:).eq.1) ! previous prescribed thickness |
---|
| 101 | vieuxH(:,:)=Hp(:,:) ! H may have been changed by calving |
---|
| 102 | elsewhere |
---|
| 103 | vieuxH(:,:)=H(:,:) |
---|
| 104 | end where |
---|
[77] | 105 | !$OMP END WORKSHARE |
---|
[4] | 106 | |
---|
| 107 | ! mk0 est le masque à ne jamais dépasser |
---|
| 108 | ! mk0=0 -> pas de glace |
---|
| 109 | ! mk0=-1 -> on garde la même epaisseur |
---|
| 110 | ! pour l'Antarctique masque mko=1 partout |
---|
| 111 | |
---|
| 112 | |
---|
| 113 | Maxdia = -1.0 |
---|
| 114 | imaxdia=0 |
---|
| 115 | jmaxdia=0 |
---|
| 116 | |
---|
[77] | 117 | ! attention limitation ! |
---|
| 118 | !$OMP WORKSHARE |
---|
[4] | 119 | uxbar(:,:)=max(uxbar(:,:),-V_limit) |
---|
| 120 | uxbar(:,:)=min(uxbar(:,:),V_limit) |
---|
| 121 | uybar(:,:)=max(uybar(:,:),-V_limit) |
---|
| 122 | uybar(:,:)=min(uybar(:,:),V_limit) |
---|
[77] | 123 | !$OMP END WORKSHARE |
---|
[4] | 124 | |
---|
| 125 | ! le pas de temps est choisi pour rendre la matrice advective diagonale dominante |
---|
[77] | 126 | !$OMP END PARALLEL |
---|
| 127 | !$OMP PARALLEL |
---|
| 128 | !$OMP DO PRIVATE(aux) |
---|
| 129 | !do i=2,nx-1 |
---|
| 130 | ! do j=2,ny-1 |
---|
| 131 | do j=2,ny-1 |
---|
| 132 | do i=2,nx-1 |
---|
[4] | 133 | aux = (abs(uxbar(i,j)) + abs(uxbar(i+1,j))) & |
---|
| 134 | + (abs(uybar(i,j)) + abs(uybar(i,j+1))) |
---|
| 135 | aux=aux*dx1*0.5 |
---|
[77] | 136 | !$OMP CRITICAL |
---|
[4] | 137 | if (aux.gt.maxdia) then |
---|
| 138 | maxdia=aux |
---|
[77] | 139 | ! imaxdia=i |
---|
| 140 | ! jmaxdia=j |
---|
[4] | 141 | endif |
---|
[77] | 142 | !$OMP END CRITICAL |
---|
[4] | 143 | end do |
---|
| 144 | end do |
---|
[77] | 145 | !$OMP END DO |
---|
| 146 | !$OMP END PARALLEL |
---|
[4] | 147 | |
---|
| 148 | if (maxdia.ge.(testdiag/dtmax)) then |
---|
| 149 | dt = testdiag/Maxdia |
---|
| 150 | dt=max(dt,dtmin) |
---|
| 151 | else |
---|
| 152 | dt = dtmax |
---|
| 153 | end if |
---|
| 154 | |
---|
| 155 | ! next_time ajuste le pas de temps pour faciliter la synchronisation |
---|
| 156 | ! avec le pas de temps long (dtt) |
---|
| 157 | |
---|
[77] | 158 | call next_time(time,dt,dtt,dtmax,dtmin,isynchro,itracebug,num_tracebug) |
---|
[4] | 159 | |
---|
| 160 | |
---|
| 161 | ! calcul diffusivite - vitesse |
---|
| 162 | !----------------------------------------------------------------- |
---|
[77] | 163 | !$OMP PARALLEL |
---|
| 164 | !$OMP WORKSHARE |
---|
[4] | 165 | Uxdiff(:,:)=diffmx(:,:)*(-sdx(:,:)) ! vitesse qui peut s'exprimer en terme diffusif |
---|
| 166 | Uydiff(:,:)=diffmy(:,:)*(-sdy(:,:)) ! dans les where qui suivent elle est limitee a uxbar |
---|
| 167 | ! pour eviter des problemes dans le cas 0< adv_frac <1 |
---|
| 168 | dmx(:,:)=diffmx(:,:) |
---|
| 169 | dmy(:,:)=diffmy(:,:) |
---|
| 170 | |
---|
| 171 | ! schema amont pour la diffusion en x |
---|
| 172 | where (Uxbar(:,:).ge.0.) |
---|
| 173 | dmx(:,:)=diffmx(:,:)*eoshift(H(:,:),shift=-1,boundary=0.,dim=1) |
---|
| 174 | uxdiff(:,:)=min(uxdiff(:,:),uxbar(:,:)) ! pour tenir compte limitation utot |
---|
| 175 | elsewhere |
---|
| 176 | dmx(:,:)=diffmx(:,:)*H(:,:) |
---|
| 177 | uxdiff(:,:)=max(uxdiff(:,:),uxbar(:,:)) ! pour tenir compte limitation utot |
---|
| 178 | end where |
---|
| 179 | |
---|
| 180 | |
---|
| 181 | |
---|
| 182 | ! schema amont pour la diffusion en y |
---|
| 183 | where (uybar(:,:).ge.0.) |
---|
| 184 | dmy(:,:)=diffmy(:,:)*eoshift(H(:,:),shift=-1,boundary=0.,dim=2) |
---|
| 185 | uydiff(:,:)=min(uydiff(:,:),uybar(:,:)) ! pour tenir compte limitation utot |
---|
| 186 | elsewhere |
---|
| 187 | dmy(:,:)=dmy(:,:)*H(:,:) |
---|
| 188 | uydiff(:,:)=max(uydiff(:,:),uybar(:,:)) ! pour tenir compte limitation utot |
---|
| 189 | end where |
---|
[77] | 190 | !$OMP END WORKSHARE |
---|
[4] | 191 | |
---|
| 192 | |
---|
| 193 | ! tests pour la répartition advection - diffusion |
---|
| 194 | !------------------------------------------------ |
---|
| 195 | |
---|
| 196 | if (adv_frac.gt.1.) then ! advection seulement |
---|
[77] | 197 | !$OMP WORKSHARE |
---|
[4] | 198 | advmx(:,:)=Uxbar(:,:) |
---|
| 199 | advmy(:,:)=Uybar(:,:) |
---|
| 200 | Dmx(:,:)=0. |
---|
| 201 | Dmy(:,:)=0. |
---|
[77] | 202 | !$OMP END WORKSHARE |
---|
[4] | 203 | else if (abs(adv_frac).lt.1.e-8) then ! diffusion seulement |
---|
[77] | 204 | !$OMP WORKSHARE |
---|
[4] | 205 | advmx(:,:)=0. |
---|
| 206 | advmy(:,:)=0. |
---|
[77] | 207 | !$OMP END WORKSHARE |
---|
[4] | 208 | ! D reste la valeur au dessus) |
---|
| 209 | |
---|
| 210 | else if ((adv_frac.ge.1.e-8).and.(adv_frac.le.1.)) then ! advection-diffusion) |
---|
| 211 | |
---|
| 212 | ! selon x -------------------------------- |
---|
| 213 | |
---|
| 214 | ! advection = adv_frac* tout ce qui n'est pas diffusion |
---|
| 215 | ! diffusion = ce qui peut être exprime en diffusion |
---|
| 216 | ! + une partie (1-adv_frac) de l'advection |
---|
[77] | 217 | !$OMP WORKSHARE |
---|
[4] | 218 | where ((abs(sdx(:,:)).gt.1.e-8).and.(Hmx(:,:).gt.2.)) ! cas general |
---|
| 219 | advmx(:,:)=(Uxbar(:,:)-Uxdiff(:,:)) ! tout ce qui n'est pas diffusion |
---|
| 220 | advmx(:,:)=advmx(:,:)*adv_frac ! la fraction adv_frac du precedent |
---|
| 221 | Dminx(:,:)=-(Uxbar(:,:)-advmx(:,:))/sdx(:,:) ! ce qui reste exprime en diffusivite |
---|
| 222 | ! a multiplier par H |
---|
| 223 | |
---|
| 224 | ! schema amont pour la diffusivite |
---|
| 225 | where (uxbar(:,:).ge.0.) |
---|
| 226 | dmx(:,:)=dminx(:,:)*eoshift(H(:,:),shift=-1,boundary=0.,dim=1) |
---|
| 227 | elsewhere |
---|
| 228 | dmx(:,:)=dminx(:,:)*H(:,:) |
---|
| 229 | end where |
---|
| 230 | |
---|
| 231 | |
---|
| 232 | elsewhere ! zones trop plates ou trop fines |
---|
| 233 | Dmx(:,:)=0. |
---|
| 234 | advmx(:,:)=Uxbar(:,:) |
---|
| 235 | end where |
---|
| 236 | |
---|
| 237 | |
---|
| 238 | ! selon y -------------------------------- |
---|
| 239 | |
---|
| 240 | ! advection = adv_frac* tout ce qui n'est pas diffusion |
---|
| 241 | ! diffusion = ce qui peut être exprime en diffusion |
---|
| 242 | ! + une partie (1-adv_frac) de l'advection |
---|
| 243 | |
---|
| 244 | where ((abs(sdy(:,:)).gt.1.e-8).and.(Hmy(:,:).gt.2.)) ! cas general |
---|
| 245 | advmy(:,:)=(Uybar(:,:)-Uydiff(:,:)) ! tout ce qui n'est pas diffusion |
---|
| 246 | advmy(:,:)=advmy(:,:)*adv_frac ! la fraction adv_frac du precedent |
---|
| 247 | Dminy(:,:)=-(Uybar(:,:)-advmy(:,:))/sdy(:,:) ! ce qui reste exprime en diffusivite |
---|
| 248 | ! a multiplier par H |
---|
| 249 | |
---|
| 250 | ! schema amont pour la diffusivite |
---|
| 251 | where (uybar(:,:).ge.0.) |
---|
| 252 | dmy(:,:)=dminy(:,:)*eoshift(H(:,:),shift=-1,boundary=0.,dim=2) |
---|
| 253 | elsewhere |
---|
| 254 | dmy(:,:)=dminy(:,:)*H(:,:) |
---|
| 255 | end where |
---|
| 256 | elsewhere ! zones trop plates ou trop fines |
---|
| 257 | Dmy(:,:)=0. |
---|
| 258 | advmy(:,:)=Uybar(:,:) |
---|
| 259 | end where |
---|
[77] | 260 | !$OMP END WORKSHARE |
---|
[4] | 261 | |
---|
| 262 | else if (adv_frac.lt.0) then ! diffusif dans zones SIA, advectif ailleurs |
---|
[77] | 263 | !$OMP WORKSHARE |
---|
[4] | 264 | zonemx(:,:)=flgzmx(:,:) |
---|
| 265 | zonemy(:,:)=flgzmy(:,:) |
---|
| 266 | |
---|
| 267 | ! selon x -------------- |
---|
| 268 | where (zonemx(:,:)) |
---|
| 269 | advmx(:,:)=Uxbar(:,:) |
---|
| 270 | Dmx(:,:)=0. |
---|
| 271 | elsewhere |
---|
| 272 | advmx(:,:)=0. |
---|
| 273 | end where |
---|
| 274 | |
---|
| 275 | ! selon y -------------- |
---|
| 276 | where (zonemy(:,:)) |
---|
| 277 | advmy(:,:)=Uybar(:,:) |
---|
| 278 | Dmy(:,:)=0. |
---|
| 279 | elsewhere |
---|
| 280 | advmy(:,:)=0. |
---|
| 281 | end where |
---|
[77] | 282 | !$OMP END WORKSHARE |
---|
[4] | 283 | end if |
---|
| 284 | |
---|
[77] | 285 | !$OMP WORKSHARE |
---|
[4] | 286 | where (flot(:,:) ) |
---|
| 287 | tabtest(:,:)=1. |
---|
| 288 | elsewhere |
---|
| 289 | tabtest(:,:)=0. |
---|
| 290 | end where |
---|
[77] | 291 | !$OMP END WORKSHARE |
---|
| 292 | !$OMP END PARALLEL |
---|
[4] | 293 | !------------------------------------------------------------------detect |
---|
| 294 | !!$ tabtest(:,:)=dmx(:,:) |
---|
| 295 | !!$ |
---|
| 296 | !!$call detect_assym(nx,ny,0,41,1,0,1,0,tabtest,itesti) |
---|
| 297 | !!$ |
---|
| 298 | !!$if (itesti.gt.0) then |
---|
| 299 | !!$ write(6,*) 'asymetrie sur dmx avant resol pour time=',time,itesti |
---|
| 300 | !!$ stop |
---|
| 301 | !!$else |
---|
| 302 | !!$ write(6,*) ' pas d asymetrie sur dmx avant resol pour time=',time,itesti |
---|
| 303 | !!$end if |
---|
| 304 | !----------------------------------------------------------- fin detect |
---|
| 305 | |
---|
| 306 | |
---|
| 307 | !------------------------------------------------------------------detect |
---|
| 308 | !!$ tabtest(:,:)=dmy(:,:) |
---|
| 309 | !!$ |
---|
| 310 | !!$call detect_assym(nx,ny,0,41,1,0,1,1,tabtest,itesti) |
---|
| 311 | !!$ |
---|
| 312 | !!$if (itesti.gt.0) then |
---|
| 313 | !!$ write(6,*) 'asymetrie sur dmy avant resol pour time=',time,itesti |
---|
| 314 | !!$ stop |
---|
| 315 | !!$else |
---|
| 316 | !!$ write(6,*) ' pas d asymetrie sur dmy avant resol pour time=',time,itesti |
---|
| 317 | !!$end if |
---|
| 318 | !----------------------------------------------------------- fin detect |
---|
| 319 | |
---|
| 320 | !------------------------------------------------------------------detect |
---|
| 321 | !!$ tabtest(:,:)=advmx(:,:) |
---|
| 322 | !!$ |
---|
| 323 | !!$call detect_assym(nx,ny,0,41,1,0,1,0,tabtest,itesti) |
---|
| 324 | !!$ |
---|
| 325 | !!$if (itesti.gt.0) then |
---|
| 326 | !!$ write(6,*) 'asymetrie sur advmx avant resol pour time=',time,itesti |
---|
| 327 | !!$ stop |
---|
| 328 | !!$else |
---|
| 329 | !!$ write(6,*) ' pas d asymetrie sur advdmx avant resol pour time=',time,itesti |
---|
| 330 | !!$end if |
---|
| 331 | !----------------------------------------------------------- fin detect |
---|
| 332 | |
---|
| 333 | !------------------------------------------------------------------detect |
---|
| 334 | !!$ tabtest(:,:)=advmy(:,:) |
---|
| 335 | !!$ |
---|
| 336 | !!$call detect_assym(nx,ny,0,41,1,0,-1,1,tabtest,itesti) |
---|
| 337 | !!$ |
---|
| 338 | !!$if (itesti.gt.0) then |
---|
| 339 | !!$ write(6,*) 'asymetrie sur advmy avant resol pour time=',time,itesti |
---|
| 340 | !!$ stop |
---|
| 341 | !!$else |
---|
| 342 | !!$ write(6,*) ' pas d asymetrie sur advmy avant resol pour time=',time,itesti |
---|
| 343 | !!$end if |
---|
| 344 | !----------------------------------------------------------- fin detect |
---|
| 345 | |
---|
| 346 | ! nouveau calcul de dt |
---|
| 347 | aux=maxval( (abs(advmx(:,:))+abs(advmy(:,:)))/dx+(abs(dmx(:,:))+abs(dmy(:,:)))/dx/dx) |
---|
| 348 | |
---|
| 349 | ! write(166,*) 'critere pour dt',time-dt,testdiag/aux,dt |
---|
| 350 | |
---|
| 351 | |
---|
| 352 | if (aux.gt.1.e-20) then |
---|
| 353 | if (testdiag/aux.lt.dt) then |
---|
| 354 | time=time-dt |
---|
| 355 | dt=testdiag/aux*4. |
---|
| 356 | if (itracebug.eq.1) call tracebug("aux tres petit,deuxieme appel a next_time dt*4") |
---|
| 357 | call next_time(time,dt,dtt,dtmax,dtmin,isynchro,itracebug,num_tracebug) |
---|
| 358 | end if |
---|
| 359 | end if |
---|
| 360 | |
---|
| 361 | |
---|
| 362 | timemax=time |
---|
| 363 | |
---|
| 364 | ! etait avant dans step |
---|
| 365 | if (time.lt.TGROUNDED) then |
---|
| 366 | MARINE=.false. |
---|
| 367 | else |
---|
| 368 | MARINE=.true. |
---|
| 369 | endif |
---|
| 370 | ! fin vient de step |
---|
| 371 | |
---|
| 372 | ! les variables dtdx et dtdx2 sont globales |
---|
| 373 | Dtdx2=Dt/(Dx**2) |
---|
| 374 | dtdx=dt/dx |
---|
| 375 | |
---|
| 376 | |
---|
| 377 | if (geoplace(1:5).eq.'mism3') then ! pas de flux sur les limites laterales |
---|
| 378 | dmy(:,2) = 0. |
---|
| 379 | dmy(:,3) = 0. |
---|
| 380 | dmy(:,ny-1) = 0. |
---|
| 381 | dmy(:,ny) = 0. |
---|
| 382 | advmy(:,2) = 0. |
---|
| 383 | advmy(:,3) = 0. |
---|
| 384 | advmy(:,ny-1) = 0. |
---|
| 385 | advmy(:,ny) = 0. |
---|
| 386 | uybar(:,2) = 0 |
---|
| 387 | uybar(:,3) = 0 |
---|
| 388 | uybar(:,ny-1) = 0 |
---|
| 389 | uybar(:,ny) = 0 |
---|
| 390 | end if |
---|
| 391 | |
---|
| 392 | |
---|
[77] | 393 | !debug_3D(:,:,45)=dmx(:,:) |
---|
| 394 | !debug_3D(:,:,46)=dmy(:,:) |
---|
| 395 | !debug_3D(:,:,47)=advmx(:,:) |
---|
| 396 | !debug_3D(:,:,48)=advmy(:,:) |
---|
[4] | 397 | |
---|
[77] | 398 | !$OMP PARALLEL |
---|
| 399 | !$OMP WORKSHARE |
---|
[4] | 400 | bilmass(:,:)=bm(:,:)-bmelt(:,:) ! surface and bottom mass balance |
---|
[77] | 401 | !$OMP END WORKSHARE |
---|
| 402 | !$OMP END PARALLEL |
---|
[4] | 403 | |
---|
| 404 | ! diverses precription de l'epaisseur en fonction de l'objectif |
---|
| 405 | !------------------------------------------------------------------- |
---|
| 406 | |
---|
| 407 | |
---|
[77] | 408 | call prescribe_fixed_points ! ceux qui sont fixes pour tout le run (bords de grille, regions) |
---|
[4] | 409 | |
---|
| 410 | |
---|
| 411 | if ((igrdline.eq.1)) then ! present grounding line |
---|
| 412 | |
---|
| 413 | ! if ((time.lt.time_gl(1)).or.(nt.lt.2)) then ATTENTION test seulement si version prescribe-H_mod.f90 |
---|
| 414 | call prescribe_present_H_gl ! prescribe present thickness at the grounding line |
---|
| 415 | ! else |
---|
| 416 | ! call prescribe_retreat |
---|
| 417 | ! endif |
---|
| 418 | end if |
---|
| 419 | |
---|
| 420 | if ((igrdline.eq.2)) then ! paleo grounding line |
---|
| 421 | call prescribe_paleo_gl_shelf |
---|
| 422 | end if |
---|
| 423 | |
---|
| 424 | if ((igrdline.eq.3)) then |
---|
| 425 | ! where (flot(:,:)) |
---|
| 426 | ! mk_gr(:,:) = 0 |
---|
| 427 | ! elsewhere |
---|
| 428 | ! mk_gr(:,:) = 1 |
---|
| 429 | ! end where |
---|
| 430 | if (itracebug.eq.1) call tracebug(" avant time_step_recul") |
---|
| 431 | call time_step_recul ! version prescribe-H-i2s_mod.f90 avec use proto_recul_mod |
---|
| 432 | if (itracebug.eq.1) call tracebug(" apres time_step_recul") |
---|
| 433 | ! call prescr_ice2sea_retreat ! version prescribe-H_mod.f90 |
---|
| 434 | endif |
---|
| 435 | |
---|
| 436 | |
---|
| 437 | !if (time.ge.t_break) then ! action sur les ice shelves |
---|
| 438 | ! call melt_ice_shelves ! ATTENTION version prescribe-H_mod.f90 |
---|
| 439 | ! call break_all_ice_shelves |
---|
| 440 | !end if |
---|
| 441 | |
---|
| 442 | |
---|
| 443 | !!$! impose une variation d'épaisseur d'apres une carte a un temps donne (a affiner plus tard pour des runs transient) |
---|
| 444 | !!$if (time.eq.t_break) then |
---|
| 445 | !!$ call lect_prescribed_deltaH |
---|
| 446 | !!$else if (time.gt.t_break) then |
---|
| 447 | !!$ call prescribe_deltaH |
---|
| 448 | !!$end if |
---|
| 449 | !!$ |
---|
| 450 | !!$if (time.eq.t_break) then ! si appele apres lect_prescribed_deltaH, cumule les effets |
---|
| 451 | !!$ call break_all_ice_shelves |
---|
| 452 | !!$end if |
---|
| 453 | |
---|
[77] | 454 | !debug_3D(:,:,87)=S(:,:)-debug_3D(:,:,88) |
---|
[4] | 455 | |
---|
[77] | 456 | !debug_3D(:,:,86)=Hp(:,:) |
---|
| 457 | !debug_3D(:,:,85)=i_Hp(:,:) |
---|
[4] | 458 | |
---|
| 459 | !!$where (i_hp(:,:).eq.1) |
---|
| 460 | !!$ vieuxh(:,:)=Hp(:,:) |
---|
| 461 | !!$endwhere |
---|
| 462 | |
---|
| 463 | |
---|
| 464 | ! Appel a la routine de resolution resol_adv_diff_2D_vect |
---|
| 465 | !------------------------------------------------------------ |
---|
| 466 | ! On ecrit les vitesses sous la forme |
---|
| 467 | ! Ux = -Dmx * sdx + advmx |
---|
| 468 | |
---|
| 469 | ! Dmx, Dmy sont les termes diffusifs de la vitesse |
---|
| 470 | ! advmx, advmy sont les termes advectifs. |
---|
| 471 | ! la repartition entre les deux depend de adv_frac. |
---|
| 472 | |
---|
| 473 | ! i_Hp est un tableau des points ou l'epaisseur est imposee a la valeur Hp |
---|
| 474 | ! bilmass est le bilan de masse (surface et fond) |
---|
| 475 | ! vieuxH est la precedente valeur de H |
---|
| 476 | ! H est le retour de la nouvelle valeur. |
---|
| 477 | |
---|
| 478 | |
---|
| 479 | call resol_adv_diff_2D_vect (Dmx,Dmy,advmx,advmy,Hp,i_Hp,bilmass,vieuxH,H) |
---|
| 480 | |
---|
[77] | 481 | ! remise à 0 des epaisseurs negatives, on garde la difference dans ablbord |
---|
| 482 | !$OMP PARALLEL |
---|
| 483 | !$OMP WORKSHARE |
---|
[4] | 484 | where (H(:,:).lt.0.) |
---|
| 485 | ablbord(:,:)=H(:,:)/dt |
---|
| 486 | elsewhere |
---|
| 487 | ablbord(:,:)=0. |
---|
| 488 | endwhere |
---|
| 489 | |
---|
| 490 | H(:,:)=max(H(:,:),0.) ! pas d'epaisseur negative |
---|
| 491 | |
---|
| 492 | |
---|
| 493 | ! calcul du masque "ice" |
---|
| 494 | where (flot(:,:)) ! points flottants, sera éventuellement réévalué dans flottab |
---|
| 495 | H(:,:)=max(H(:,:),1.) ! dans la partie marine l'épaisseur mini est 1 m |
---|
| 496 | where(H(:,:).gt.1.) |
---|
| 497 | ice(:,:)=1 |
---|
| 498 | elsewhere |
---|
| 499 | ice(:,:)=0 |
---|
| 500 | endwhere |
---|
| 501 | elsewhere ! points posés |
---|
| 502 | where(H(:,:).gt.0.) |
---|
| 503 | ice(:,:)=1 |
---|
| 504 | elsewhere |
---|
| 505 | ice(:,:)=0 |
---|
| 506 | endwhere |
---|
| 507 | endwhere |
---|
| 508 | |
---|
| 509 | ! eventuellement retirer apres spinup ou avoir un cas serac |
---|
| 510 | Hdot(:,:)=(H(:,:)-vieuxH(:,:))/dt |
---|
[77] | 511 | !$OMP END WORKSHARE |
---|
[4] | 512 | |
---|
| 513 | if (igrdline.ne.3) then |
---|
[77] | 514 | !$OMP WORKSHARE |
---|
[4] | 515 | Hdot(:,:)=min(Hdot(:,:),10.) |
---|
| 516 | Hdot(:,:)=max(Hdot(:,:),-10.) |
---|
[77] | 517 | !$OMP END WORKSHARE |
---|
[4] | 518 | endif |
---|
| 519 | |
---|
[77] | 520 | !$OMP WORKSHARE |
---|
[4] | 521 | where (i_hp(:,:).ne.1) |
---|
| 522 | H(:,:)=vieuxH(:,:)+Hdot(:,:)*dt |
---|
| 523 | end where |
---|
| 524 | |
---|
| 525 | |
---|
| 526 | |
---|
| 527 | ! si mk0=-1, on garde l'epaisseur precedente |
---|
| 528 | ! en attendant un masque plus propre |
---|
| 529 | |
---|
| 530 | where(mk0(:,:).eq.-1) |
---|
| 531 | H(:,:)=vieuxH(:,:) |
---|
| 532 | Hdot(:,:)=0. |
---|
| 533 | end where |
---|
[77] | 534 | !$OMP END WORKSHARE |
---|
| 535 | !$OMP END PARALLEL |
---|
[4] | 536 | |
---|
| 537 | !calul de l'ablation sur les bords (pourrait n'être appelé que pour les sorties courtes) |
---|
| 538 | if (isynchro.eq.1) call ablation_bord |
---|
| 539 | |
---|
| 540 | |
---|
| 541 | if (itracebug.eq.1) call tracebug(' Fin routine icethick') !, maxval(H) ',maxval(H(:,:)) |
---|
| 542 | |
---|
| 543 | |
---|
| 544 | end subroutine icethick3 |
---|
| 545 | |
---|
| 546 | end module equat_adv_diff_2D_vect |
---|