Changeset 237 for trunk/SOURCES/Eurasie40_files
- Timestamp:
- 01/09/19 17:09:26 (5 years ago)
- Location:
- trunk/SOURCES/Eurasie40_files
- Files:
-
- 4 edited
Legend:
- Unmodified
- Added
- Removed
-
trunk/SOURCES/Eurasie40_files/bmelt-eurasie-depth-lake_mod.f90
r4 r237 85 85 do i=1,nx 86 86 87 if ((sealevel -BSOC(i,j)).lt.450) then87 if ((sealevel_2d(i,j)-BSOC(i,j)).lt.450) then 88 88 bmshelf(i,j)=2.00 !3.00! 2.55! 0.45! 0.65 89 89 bmgrz(i,j)=2.00 !3.00! 2.55! 0.45 ! 0.65 … … 127 127 128 128 129 if (s lv(i,j).gt.0) then129 if (sealevel_2d(i,j).gt.0) then 130 130 bmshelf(i,j)=0.05!5.0 131 131 bmgrz(i,j)=0.05!5.0 … … 135 135 136 136 137 if (s lv(i,j).gt.sealevel) then137 if (sealevel_2d(i,j).gt.sealevel) then 138 138 bmelt(i,j)=bmshelf(i,j) 139 139 if (fbm(i,j)) bmelt(i,j)=bmgrz(i,j) … … 166 166 ! en fonction du nombre de points flottants 167 167 168 if (s lv(i,j).gt.sealevel) then168 if (sealevel_2d(i,j).gt.sealevel) then 169 169 bmelt(i,j)= ngr/4.*bmgrz(i,j)+(1.-ngr/4.)*bmelt(i,j) 170 170 else -
trunk/SOURCES/Eurasie40_files/climat-forcage-eurasie_mod-0.4.f90
r4 r237 493 493 nom_table='lvflo' 494 494 call printtable_r(levelflot,nom_table) 495 495 496 sealevel_2d(:,:) = sealevel 497 496 498 ! coefbmshelf(:,:)=(1.+delTatime(:,:)/7.) 497 499 !coefbmshelf=1. ! modif pour test de sensibilite (tof 20 avril 01) -
trunk/SOURCES/Eurasie40_files/lakes-prescribed_mod-0.1.f90
r4 r237 54 54 subroutine lake_to_slv 55 55 56 ! inclut les lacs dans le tableau s lvdes niveaux de flottaison56 ! inclut les lacs dans le tableau sealevel_2d des niveaux de flottaison 57 57 58 s lv(:,:)=sealevel58 sealevel_2d(:,:)=sealevel 59 59 60 60 ! VOIR COMMENT ACTIVER-DESACTIVER LES LACS EN FONCTION DU TEMPS … … 64 64 do i_lake=1,nb_lakes 65 65 if (lake(i,j).eq.i_lake) then 66 s lv(i,j)=max(level_lakes(i_lake),Bsoc(i,j))66 sealevel_2d(i,j)=max(level_lakes(i_lake),Bsoc(i,j)) 67 67 end if 68 68 end do -
trunk/SOURCES/Eurasie40_files/output_eurasie40_mod-0.1.f90
r4 r237 85 85 86 86 ! calcul de la hauteur au dessus de la flottaison 87 if (sealevel -B(i,j).le.0.) then ! socle au dessus du niveau des mers87 if (sealevel_2d(i,j)-B(i,j).le.0.) then ! socle au dessus du niveau des mers 88 88 volf=volf+h(i,j) 89 89 else 90 volf=volf+h(i,j)-row/ro*(sealevel -b(i,j))90 volf=volf+h(i,j)-row/ro*(sealevel_2d(i,j)-b(i,j)) 91 91 endif 92 92
Note: See TracChangeset
for help on using the changeset viewer.