#ifdef ISOVERIF ! $Id: $ MODULE isotopes_verif_mod !use isotopes_mod, ONLY: !#ifdef ISOTRAC !use isotrac_mod, ONLY: !#endif implicit none save ! variables de vérifications real errmax ! erreur maximale en absolu. parameter (errmax=1e-8) real errmax_sol ! erreur maximale en absolu. parameter (errmax_sol=5e-7) real errmaxrel ! erreur maximale en relatif autorisée ! parameter (errmaxrel=1e10) parameter (errmaxrel=1e-3) real borne ! valeur maximale que n'importe quelle variable peut ! atteindre, en val abs; utile pour vérif que pas NAN parameter (borne=1e19) real, save :: deltalim ! deltalim est le maximum de deltaD qu'on ! autorise dans la vapeur, au-dela, en suppose que c'est abérrant. ! dans le liquide, on autorise deltalim*faccond. !$OMP THREADPRIVATE(deltalim) ! parameter (deltalim=1e10) ! parameter (deltalim=300.0) ! maintenant défini dans iso.def real, save :: deltalimtrac ! max de deltaD dans les traceurs, défini dans iso.def !$OMP THREADPRIVATE(deltalimtrac) real, save :: deltalim_snow ! deltalim est le maximum de deltaD qu'on ! autorise dans la neige, au-dela, en suppose que c'est abérrant. !$OMP THREADPRIVATE(deltalim_snow) ! parameter (deltalim_snow=500.0) ! initalisé dans iso_init real, save :: deltaDmin !$OMP THREADPRIVATE(deltaDmin) ! parameter (deltaDmin=-900.0) ! maintentant, défini dans iso.def real, save :: o17excess_bas,o17excess_haut ! bornes inf et sup de l'O17-excess ! parameter(o17excess_bas=-200.0,o17excess_haut=120) !$OMP THREADPRIVATE(o17excess_bas,o17excess_haut) integer nlevmaxO17 !$OMP THREADPRIVATE(nlevmaxO17) logical, save :: O17_verif !$OMP THREADPRIVATE(O17_verif) ! parameter (O17_verif=.true.) real, save :: dexcess_min,dexcess_max !$OMP THREADPRIVATE(dexcess_min,dexcess_max) real faccond ! dans le liquide, on autorise R(deltalim)*faccond. parameter (faccond=1.1) ! logical bidouille_anti_divergence ! si true, alors on fait un ! ! rappel régulier des xt4 vers les q pour ! !éviter accumulation d'érreurs et divergences ! parameter (bidouille_anti_divergence=.true.) ! parameter (bidouille_anti_divergence=.false.) ! bidouille_anti_divergence a été déplacé dans wateriso2.h et est lue à l'execution real deltaDfaible ! deltaD en dessous duquel la vapeur est tellement faible !que on peut accepter que la remise à l'équilibre du sol avec ! cette vapeur donne des deltaDevap aberrants. parameter (deltaDfaible=-300) real deltaDfaible_lax ! deltaD en dessous duquel la vapeur est tellement faible !que on peut accepter que la remise à l'équilibre du sol avec ! cette vapeur donne des deltaDevap aberrants. parameter (deltaDfaible_lax=-180) real faible_evap ! faible évaporation: on est plus laxiste !pour les deltaD aberrant dans le cas de l'évap venant d'orchidee parameter (faible_evap=3.0) real Tmin_verif parameter (Tmin_verif=5.0) ! en K real Tmax_verif parameter (Tmax_verif=1000.0) ! en K ! subroutines de vérifications génériques, à ne pas modifier CONTAINS SUBROUTINE iso_verif_init() use ioipsl_getin_p_mod, ONLY : getin_p !USE infotrac_phy, ONLY: use_iso use isotopes_mod, ONLY: iso_O17, iso_O18, iso_HDO implicit none write(*,*) 'iso_verif_init 99: entree' deltalim=300.0 deltalim_snow=500.0 deltaDmin=-900.0 deltalimtrac=2000.0 O17_verif=.false. o17excess_bas=-200.0 o17excess_haut=120.0 dexcess_min=-100.0 dexcess_max=1000.0 call getin_p('deltalim',deltalim) deltalim_snow=deltalim+200.0 call getin_p('deltaDmin',deltaDmin) call getin_p('deltalimtrac',deltalimtrac) write(*,*) 'iso_O17,iso_O18,iso_HDO=',iso_O17,iso_O18,iso_HDO if ((iso_O17.gt.0).and.(iso_O18.gt.0)) then call getin_p('O17_verif',O17_verif) call getin_p('o17excess_bas',o17excess_bas) call getin_p('o17excess_haut',o17excess_haut) call getin_p('nlevmaxO17',nlevmaxO17) endif if ((iso_HDO.gt.0).and.(iso_O18.gt.0)) then call getin_p('dexcess_min',dexcess_min) call getin_p('dexcess_max',dexcess_max) endif write(*,*) 'deltalim=',deltalim write(*,*) 'deltaDmin=',deltaDmin write(*,*) 'deltalimtrac=',deltalimtrac write(*,*) 'O17_verif=',O17_verif write(*,*) 'o17excess_bas=',o17excess_bas write(*,*) 'o17excess_haut=',o17excess_haut write(*,*) 'dexcess_min=',dexcess_min write(*,*) 'dexcess_max=',dexcess_max end SUBROUTINE iso_verif_init subroutine iso_verif_egalite(a,b,err_msg) implicit none ! compare a et b. Si pas egal à errmax près, on affiche message ! d'erreur et stoppe ! input: real a, b character*(*) err_msg ! message d''erreur à afficher ! local !integer iso_verif_egalite_choix_nostop if (iso_verif_egalite_choix_nostop(a,b,err_msg, & & errmax,errmaxrel).eq.1) then stop endif #ifdef ISOVERIF #else write(*,*) 'err_msg=',err_msg,': pourquoi verif?' stop #endif return end subroutine iso_verif_egalite !***************** function iso_verif_egalite_nostop(a,b,err_msg) implicit none ! compare a et b. Si pas egal à errmax près, on affiche message ! d'erreur et stoppe ! input: real a, b character*(*) err_msg ! message d''erreur à afficher !ouptut integer iso_verif_egalite_nostop ! local !integer iso_verif_egalite_choix_nostop iso_verif_egalite_nostop=iso_verif_egalite_choix_nostop & & (a,b,err_msg,errmax,errmaxrel) #ifdef ISOVERIF #else write(*,*) 'err_msg=',err_msg,': pourquoi verif?' stop #endif return end function iso_verif_egalite_nostop !************************************ subroutine iso_verif_aberrant(R,err_msg) !USE infotrac_phy, ONLY: use_iso use isotopes_mod, ONLY: ridicule, iso_HDO implicit none ! si le rapprot iso R est plus grand que deltaDlim, on affiche message ! d'erreur et stoppe ! input: real R character*(*) err_msg ! message d''erreur à afficher !real deltaD !integer iso_verif_aberrant_choix_nostop if (iso_verif_aberrant_choix_nostop(R,1.0,ridicule, & & deltalim,err_msg).eq.1) then stop endif #ifdef ISOVERIF if (.not.(iso_HDO.gt.0)) then write(*,*) 'iso86: err_msg=',err_msg,': pourquoi verif?' stop endif #else write(*,*) 'iso 90: err_msg=',err_msg,': pourquoi verif?' stop #endif end subroutine iso_verif_aberrant subroutine iso_verif_aberrant_encadre(R,err_msg) !use infotrac_phy, ONLY: use_iso use isotopes_mod, ONLY: ridicule, iso_HDO implicit none ! si le rapprot iso R est plus grand que deltaDlim, on affiche message ! d'erreur et stoppe ! input: real R character*(*) err_msg ! message d''erreur à afficher !real deltaD !integer iso_verif_aberrant_enc_choix_nostop if (iso_verif_aberrant_enc_choix_nostop(R,1.0,ridicule, & & deltalim,err_msg).eq.1) then write(*,*) 'deltaD=',deltaD(R) call abort_physic('isotopes_verif_mod > iso_verif_aberrant_encadre',err_msg,1) !stop endif #ifdef ISOVERIF if (.not.(iso_HDO.gt.0)) then write(*,*) 'iso86: err_msg=',err_msg,': pourquoi verif?' stop endif #else write(*,*) 'iso 90: err_msg=',err_msg,': pourquoi verif?' stop #endif end subroutine iso_verif_aberrant_encadre !************************************ subroutine iso_verif_aberrant_choix(xt,q,qmin,deltaDmax,err_msg) !use infotrac_phy, ONLY: use_iso use isotopes_mod, ONLY: iso_HDO implicit none ! si le rapprot iso R est plus grand que deltaDlim, on affiche message ! d'erreur et stoppe ! input: real xt,q,qmin,deltaDmax character*(*) err_msg ! message d''erreur à afficher !real deltaD ! locals !integer iso_verif_aberrant_choix_nostop if (iso_verif_aberrant_choix_nostop(xt,q,qmin, & & deltaDmax,err_msg).eq.1) then stop endif #ifdef ISOVERIF if (.not.(iso_HDO.gt.0)) then write(*,*) 'iso122: err_msg=',err_msg,': pourquoi verif?' stop endif #else write(*,*) 'iso126: err_msg=',err_msg,': pourquoi verif?' stop #endif end subroutine iso_verif_aberrant_choix !************************************ function iso_verif_aberrant_nostop(R,err_msg) !use infotrac_phy, ONLY: use_iso use isotopes_mod, ONLY: ridicule,iso_HDO implicit none ! si le rapprot iso R est plus grand que deltaDlim, on affiche message ! d'erreur et stoppe ! input: real R character*(*) err_msg ! message d''erreur à afficher integer iso_verif_aberrant_nostop ! output: 1 si erreur, 0 sinon !real deltaD ! locals !integer iso_verif_aberrant_choix_nostop iso_verif_aberrant_nostop=iso_verif_aberrant_choix_nostop & & (R,1.0,ridicule,deltalim,err_msg) #ifdef ISOVERIF if (.not.(iso_HDO.gt.0)) then write(*,*) 'iso156: err_msg=',err_msg,': pourquoi verif?' stop endif #else write(*,*) 'iso160: err_msg=',err_msg,': pourquoi verif?' stop #endif return end function iso_verif_aberrant_nostop function iso_verif_aberrant_enc_nostop(R,err_msg) !use infotrac_phy, ONLY: use_iso use isotopes_mod, ONLY: ridicule,iso_HDO implicit none ! si le rapprot iso R est plus grand que deltaDlim, on affiche message ! d'erreur et stoppe ! input: real R character*(*) err_msg ! message d''erreur à afficher integer iso_verif_aberrant_enc_nostop ! output: 1 si erreur, 0 sinon !real deltaD ! locals !integer iso_verif_aberrant_enc_choix_nostop iso_verif_aberrant_enc_nostop= & & iso_verif_aberrant_enc_choix_nostop & & (R,1.0,ridicule,deltalim,err_msg) #ifdef ISOVERIF if (.not.(iso_HDO.gt.0)) then write(*,*) 'iso156: err_msg=',err_msg,': pourquoi verif?' stop endif #else write(*,*) 'iso160: err_msg=',err_msg,': pourquoi verif?' stop #endif return end function iso_verif_aberrant_enc_nostop !************************************ function iso_verif_aberrant_choix_nostop(xt,q, & & qmin,deltaDmax,err_msg) !use infotrac_phy, ONLY: use_iso use isotopes_mod, ONLY: iso_HDO implicit none ! si le rapprot iso R est plus grand que deltaDlim, on affiche message ! d'erreur et stoppe ! input: real xt,q,qmin,deltaDmax character*(*) err_msg ! message d''erreur à afficher ! output integer iso_verif_aberrant_choix_nostop ! locals !real deltaD !integer iso_verif_noNaN_nostop iso_verif_aberrant_choix_nostop=0 #ifdef ISOVERIF if (iso_verif_noNaN_nostop(q,err_msg).eq.1) then write(*,*) 'q=',q iso_verif_aberrant_choix_nostop=1 endif if (iso_verif_noNaN_nostop(xt,err_msg).eq.1) then write(*,*) 'xt=',xt iso_verif_aberrant_choix_nostop=1 endif #endif if (q.gt.qmin) then if ((deltaD(xt/q).gt.deltaDmax).or. & & (deltaD(xt/q).lt.-borne).or. & & (xt.lt.-borne).or. & & (xt.gt.borne)) then write(*,*) 'erreur detectee par '// & & 'iso_verif_aberrant_choix_nostop:' write(*,*) err_msg write(*,*) 'q,deltaD=',q,deltaD(xt/q) write(*,*) 'deltaDmax=',deltaDmax iso_verif_aberrant_choix_nostop=1 if (abs(xt-q)/q.lt.errmax) then write(*,*) 'attention, n''a-t-on pas confondu'// & & ' iso_HDO et iso_eau dans l''appel à la verif?' endif endif endif #ifdef ISOVERIF if (.not.(iso_HDO.gt.0)) then write(*,*) 'iso205: err_msg=',err_msg,': pourquoi verif?' stop endif #else write(*,*) 'iso 209: err_msg=',err_msg,': pourquoi verif?' stop #endif return end function iso_verif_aberrant_choix_nostop function iso_verif_aberrant_enc_choix_nostop(xt,q, & & qmin,deltaDmax,err_msg) !use infotrac_phy, ONLY: use_iso use isotopes_mod, ONLY: iso_HDO implicit none ! si le rapprot iso R est plus grand que deltaDlim, on affiche message ! d'erreur et stoppe ! input: real xt,q,qmin,deltaDmax character*(*) err_msg ! message d''erreur à afficher ! output integer iso_verif_aberrant_enc_choix_nostop ! locals !real deltaD iso_verif_aberrant_enc_choix_nostop=0 if (q.gt.qmin) then if ((deltaD(xt/q).gt.deltaDmax).or. & & (deltaD(xt/q).lt.deltaDmin)) then write(*,*) 'erreur detectee par '// & & 'iso_verif_aberrant_choix_nostop:' write(*,*) err_msg write(*,*) 'q,deltaD=',q,deltaD(xt/q) iso_verif_aberrant_enc_choix_nostop=1 if (abs(xt-q)/q.lt.errmax) then write(*,*) 'attention, n''a-t-on pas confondu'// & & ' iso_HDO et iso_eau dans l''appel à la verif?' endif endif endif #ifdef ISOVERIF if (.not.(iso_HDO.gt.0)) then write(*,*) 'iso205: err_msg=',err_msg,': pourquoi verif?' stop endif #else write(*,*) 'iso 209: err_msg=',err_msg,': pourquoi verif?' stop #endif return end function iso_verif_aberrant_enc_choix_nostop !******************* subroutine iso_verif_aberrant_o17(R17,R18,err_msg) implicit none ! si l'O17-excess est aberrant, on affiche un message ! d'erreur et stoppe ! input: real R17,R18 character*(*) err_msg ! message d''erreur à afficher !real o17excess ! locals !integer iso_verif_aberrant_o17_nostop ! write(*,*) 'O17_verif=',O17_verif if (O17_verif) then if (iso_verif_aberrant_o17_nostop(R17,R18,err_msg) & & .eq.1) then stop endif endif !if (O17_verif) then #ifdef ISOVERIF #else write(*,*) 'err_msg=',err_msg,': pourquoi verif?' stop #endif return end subroutine iso_verif_aberrant_o17 !******************* function iso_verif_aberrant_o17_nostop(R17,R18,err_msg) USE isotopes_mod, ONLY: tnat,iso_O17,iso_O18 implicit none ! si l'O17-excess est aberrant, on affiche un message ! d'erreur et renvoit 1 ! input: real R17,R18 character*(*) err_msg ! message d''erreur à afficher !local !real o17excess ! output integer iso_verif_aberrant_o17_nostop if (O17_verif) then iso_verif_aberrant_o17_nostop=0 if ((o17excess(R17,R18).gt.o17excess_haut).or. & & (o17excess(R17,R18).lt.o17excess_bas)) then write(*,*) 'erreur detectee par iso_verif_aberrant_O17:' write(*,*) err_msg write(*,*) 'o17excess=',o17excess(R17,R18) write(*,*) 'deltaO17=',(R17/tnat(iso_o17)-1.0)*1000.0 write(*,*) 'deltaO18=',(R18/tnat(iso_O18)-1.0)*1000.0 ! attention, vérifier que la ligne suivante est bien activée iso_verif_aberrant_o17_nostop=1 endif endif !if (O17_verif) then #ifdef ISOVERIF #else write(*,*) 'err_msg=',err_msg,': pourquoi verif?' stop #endif return end function iso_verif_aberrant_o17_nostop !************************************ subroutine iso_verif_noNaN(x,err_msg) implicit none ! si x est NaN, on affiche message ! d'erreur et stoppe ! input: real x character*(*) err_msg ! message d''erreur à afficher ! locals !integer iso_verif_noNAN_nostop if (iso_verif_noNAN_nostop(x,err_msg).eq.1) then stop endif #ifdef ISOVERIF #else write(*,*) 'err_msg iso443=',err_msg,': pourquoi verif?' stop #endif end subroutine iso_verif_noNaN !************************************ function iso_verif_noNaN_nostop(x,err_msg) implicit none ! si x est NaN, on affiche message ! d'erreur et return 1 si erreur ! input: real x character*(*) err_msg ! message d''erreur à afficher ! output integer iso_verif_noNAN_nostop if ((x.gt.-borne).and.(x.lt.borne)) then iso_verif_noNAN_nostop=0 else write(*,*) 'erreur detectee par iso_verif_nonNAN:' write(*,*) err_msg write(*,*) 'x=',x iso_verif_noNAN_nostop=1 endif #ifdef ISOVERIF #else write(*,*) 'err_msg iso 482=',err_msg,': pourquoi verif?' stop #endif return end function iso_verif_noNaN_nostop subroutine iso_verif_noNaN_vect2D(x,err_msg,ni,n,m) implicit none ! si x est NaN, on affiche message ! d'erreur et return 1 si erreur ! input: integer n,m,ni real x(ni,n,m) character*(*) err_msg ! message d''erreur à afficher ! output ! locals integer i,j,ixt do i=1,n do j=1,m do ixt=1,ni if ((x(ixt,i,j).gt.-borne).and. & & (x(ixt,i,j).lt.borne)) then else !if ((x(ixt,i,j).gt.-borne).and. write(*,*) 'erreur detectee par iso_verif_nonNAN:' write(*,*) err_msg write(*,*) 'x,ixt,i,j=',x(ixt,i,j),ixt,i,j stop endif !if ((x(ixt,i,j).gt.-borne).and. enddo !do ixt=1,ni enddo !do j=1,m enddo !do i=1,n #ifdef ISOVERIF #else write(*,*) 'err_msg iso525=',err_msg,': pourquoi verif?' stop #endif end subroutine iso_verif_noNaN_vect2D subroutine iso_verif_noNaN_vect1D(x,err_msg,ni,n) implicit none ! si x est NaN, on affiche message ! d'erreur et return 1 si erreur ! input: integer n,ni real x(ni,n) character*(*) err_msg ! message d''erreur à afficher ! output ! locals integer i,ixt do i=1,n do ixt=1,ni if ((x(ixt,i).gt.-borne).and. & & (x(ixt,i).lt.borne)) then else !if ((x(ixt,i,j).gt.-borne).and. write(*,*) 'erreur detectee par iso_verif_nonNAN:' write(*,*) err_msg write(*,*) 'x,ixt,i=',x(ixt,i),ixt,i stop endif !if ((x(ixt,i,j).gt.-borne).and. enddo !do ixt=1,ni enddo !do i=1,n #ifdef ISOVERIF #else write(*,*) 'err_msg iso525=',err_msg,': pourquoi verif?' stop #endif end subroutine iso_verif_noNaN_vect1D !************************ subroutine iso_verif_egalite_choix(a,b,err_msg,erabs,errel) implicit none ! compare a et b. Si pas egal, on affiche message ! d'erreur et stoppe ! pour egalite, on verifie erreur absolue et arreur relative ! input: real a, b real erabs,errel !erreur absolue et relative character*(*) err_msg ! message d''erreur à afficher ! locals !integer iso_verif_egalite_choix_nostop if (iso_verif_egalite_choix_nostop( & & a,b,err_msg,erabs,errel).eq.1) then stop endif #ifdef ISOVERIF #else write(*,*) 'err_msg=',err_msg,': pourquoi verif?' stop #endif end subroutine iso_verif_egalite_choix !************************ function iso_verif_egalite_choix_nostop & & (a,b,err_msg,erabs,errel) implicit none ! compare a et b. Si pas egal, on affiche message ! d'erreur et stoppe ! pour egalite, on verifie erreur absolue et arreur relative ! input: real a, b real erabs,errel !erreur absolue et relative character*(*) err_msg ! message d''erreur à afficher ! output integer iso_verif_egalite_choix_nostop ! locals !integer iso_verif_noNaN_nostop iso_verif_egalite_choix_nostop=0 #ifdef ISOVERIF if (iso_verif_noNaN_nostop(a,err_msg).eq.1) then write(*,*) 'a=',a iso_verif_egalite_choix_nostop=1 endif if (iso_verif_noNaN_nostop(b,err_msg).eq.1) then write(*,*) 'b=',b iso_verif_egalite_choix_nostop=1 endif #endif if (abs(a-b).gt.erabs) then if (abs((a-b)/max(max(abs(b),abs(a)),1e-18)) & & .gt.errel) then write(*,*) 'erreur detectee par iso_verif_egalite:' write(*,*) err_msg write(*,*) 'a=',a write(*,*) 'b=',b write(*,*) 'erabs,errel=',erabs,errel iso_verif_egalite_choix_nostop=1 ! stop endif endif #ifdef ISOVERIF #else write(*,*) 'err_msg=',err_msg,': pourquoi verif?' stop #endif return end function iso_verif_egalite_choix_nostop !****************** subroutine iso_verif_positif(x,err_msg) use isotopes_mod, ONLY: ridicule implicit none ! si x<0, on plante. ! si très limite, on le met à 0. ! input: real x character*(*) err_msg ! message d''erreur à afficher ! locals !integer iso_verif_positif_choix_nostop if (iso_verif_positif_choix_nostop(x,ridicule,err_msg) & & .eq.1) then stop endif #ifdef ISOVERIF #else write(*,*) 'iso_verif 637: err_msg=',err_msg, & & ': pourquoi verif?' stop #endif end subroutine iso_verif_positif !****************** subroutine iso_verif_positif_vect(x,n,err_msg) use isotopes_mod, ONLY: ridicule implicit none ! si x<0, on plante. ! input: integer n real x(n) character*(*) err_msg ! message d''erreur à afficher ! locals integer i integer iso_verif_positif_nostop integer ifaux iso_verif_positif_nostop=0 do i=1,n if (x(i).lt.-ridicule) then iso_verif_positif_nostop=1 ifaux=i endif enddo if (iso_verif_positif_nostop.eq.1) then write(*,*) 'erreur detectee par iso_verif_positif_vect:' write(*,*) err_msg write(*,*) 'i,x=',ifaux,x(ifaux) stop endif end subroutine iso_verif_positif_vect subroutine iso_verif_positif_choix_vect(x,n,ridic,err_msg) implicit none ! si x<0, on plante. ! input: integer n real x(n) real ridic character*(*) err_msg ! message d''erreur à afficher ! locals integer i integer iso_verif_positif_nostop integer ifaux iso_verif_positif_nostop=0 do i=1,n if (x(i).lt.-ridic) then iso_verif_positif_nostop=1 ifaux=i endif enddo if (iso_verif_positif_nostop.eq.1) then write(*,*) 'erreur detectee par iso_verif_positif_choix_vect:' write(*,*) err_msg write(*,*) 'i,x=',ifaux,x(ifaux) stop endif end subroutine iso_verif_positif_choix_vect !****************** subroutine iso_verif_positif_strict(x,err_msg) implicit none ! si x<0, on plante. ! si très limite, on le met à 0. ! input: real x character*(*) err_msg ! message d''erreur à afficher ! locals !integer iso_verif_positif_strict_nostop if (iso_verif_positif_strict_nostop(x,err_msg) & & .eq.1) then stop endif end subroutine iso_verif_positif_strict !****************** function iso_verif_positif_strict_nostop(x,err_msg) implicit none ! si x<0, on plante. ! si très limite, on le met à 0. ! input: real x character*(*) err_msg ! message d''erreur à afficher* ! output integer iso_verif_positif_strict_nostop if (x.gt.0.0) then iso_verif_positif_strict_nostop=0 else write(*,*) 'erreur detectee par iso_verif_positif_strict:' write(*,*) err_msg write(*,*) 'x=',x iso_verif_positif_strict_nostop=1 ! stop endif return end function iso_verif_positif_strict_nostop !****************** subroutine iso_verif_positif_choix(x,ridic,err_msg) implicit none ! si x<0, on plante. ! si très limite, on le met à 0. ! input: real x real ridic character*(*) err_msg ! message d''erreur à afficher ! locals !integer iso_verif_positif_choix_nostop if (iso_verif_positif_choix_nostop(x,ridic,err_msg) & & .eq.1) then stop endif #ifdef ISOVERIF #else write(*,*) 'iso_verif 801: err_msg=',err_msg, & & ': pourquoi verif?' stop #endif end subroutine iso_verif_positif_choix !****************** function iso_verif_positif_nostop(x,err_msg) use isotopes_mod, ONLY: ridicule implicit none ! si x<0, on plante. ! si très limite, on le met à 0. ! input: real x character*(*) err_msg ! message d''erreur à afficher ! output integer iso_verif_positif_nostop ! locals !integer iso_verif_positif_choix_nostop iso_verif_positif_nostop=iso_verif_positif_choix_nostop & & (x,ridicule,err_msg) #ifdef ISOVERIF #else write(*,*) 'iso_verif 837: err_msg=',err_msg, & & ': pourquoi verif?' stop #endif return end function iso_verif_positif_nostop !****************** function iso_verif_positif_choix_nostop(x,ridic,err_msg) implicit none ! si x<0, on plante. ! si très limite, on le met à 0. ! input: real x real ridic character*(*) err_msg ! message d''erreur à afficher ! output integer iso_verif_positif_choix_nostop if (x.lt.-ridic) then write(*,*) 'erreur detectee par iso_verif_positif_nostop:' write(*,*) err_msg write(*,*) 'x=',x iso_verif_positif_choix_nostop=1 else ! x=max(x,0.0) iso_verif_positif_choix_nostop=0 endif #ifdef ISOVERIF #else write(*,*) 'iso_verif 877: err_msg=',err_msg, & & ': pourquoi verif?' stop #endif return end function iso_verif_positif_choix_nostop !************** subroutine iso_verif_O18_aberrant(Rd,Ro,err_msg) implicit none ! vérifie que: ! 1) deltaD et deltaO18 dans bonne gamme ! 2) dexcess dans bonne gamme ! input: real Rd,Ro character*(*) err_msg ! message d''erreur à afficher ! local !integer iso_verif_O18_aberrant_nostop if (iso_verif_O18_aberrant_nostop(Rd,Ro,err_msg).eq.1) then stop endif end subroutine iso_verif_O18_aberrant function iso_verif_O18_aberrant_nostop(Rd,Ro,err_msg) USE isotopes_mod, ONLY: tnat, iso_HDO, iso_O18 implicit none ! vérifie que: ! 1) deltaD et deltaO18 dans bonne gamme ! 2) dexcess dans bonne gamme ! input: real Rd,Ro character*(*) err_msg ! message d''erreur à afficher ! outputs integer iso_verif_O18_aberrant_nostop !locals real deltaD,deltao,dexcess deltaD=(Rd/tnat(iso_hdo)-1)*1000 deltao=(Ro/tnat(iso_O18)-1)*1000 dexcess=deltaD-8*deltao iso_verif_O18_aberrant_nostop=0 if ((deltaD.lt.deltaDmin).or.(deltao.lt.deltaDmin/2.0).or. & & (deltaD.gt.deltalim).or.(deltao.gt.deltalim/8.0).or. & & ((deltaD.gt.-500.0).and.((dexcess.lt.dexcess_min) & & .or.(dexcess.gt.dexcess_max)))) then write(*,*) 'erreur detectee par iso_verif_O18_aberrant:' write(*,*) err_msg write(*,*) 'delta180=',deltao write(*,*) 'deltaD=',deltaD write(*,*) 'Dexcess=',dexcess ! stop iso_verif_O18_aberrant_nostop=1 endif #ifdef ISOVERIF #else write(*,*) 'err_msg=',err_msg,': pourquoi verif?' stop #endif return end function iso_verif_O18_aberrant_nostop ! ********** function deltaD(R) !use infotrac_phy, ONLY: use_iso USE isotopes_mod, ONLY: tnat,iso_HDO implicit none real R,deltaD if (iso_HDO.gt.0) then deltaD=(R/tnat(iso_HDO)-1)*1000.0 else write(*,*) 'iso_verif_egalite>deltaD 260: iso_HDO.gt.0=', & & iso_HDO.gt.0 endif return end function deltaD ! ********** function deltaO(R) !use infotrac_phy, ONLY: use_iso USE isotopes_mod, ONLY: tnat,iso_O18 implicit none real R,deltaO if (iso_O18.gt.0) then deltaO=(R/tnat(iso_O18)-1)*1000.0 else write(*,*) 'iso_verif_egalite>deltaO18 260: iso_O18.gt.0=', & & iso_O18.gt.0 endif return end function deltaO ! ********** function dexcess(RD,RO) !use infotrac_phy, ONLY: use_iso USE isotopes_mod, ONLY: tnat,iso_O18,iso_HDO implicit none real RD,RO,deltaD,deltaO,dexcess if ((iso_O18.gt.0).and.(iso_HDO.gt.0)) then deltaD=(RD/tnat(iso_HDO)-1)*1000.0 deltaO=(RO/tnat(iso_O18)-1)*1000.0 dexcess=deltaD-8*deltaO else write(*,*) 'iso_verif_egalite 1109: iso_O18,iso_HDO=',iso_O18,iso_HDO endif return end function dexcess ! ********** function delta_all(R,ixt) USE isotopes_mod, ONLY: tnat implicit none real R,delta_all integer ixt delta_all=(R/tnat(ixt)-1)*1000.0 return end function delta_all ! ********** function delta_to_R(delta,ixt) USE isotopes_mod, ONLY: tnat implicit none real delta,delta_to_R integer ixt delta_to_R=(delta/1000.0+1.0)*tnat(ixt) return end function delta_to_R ! ********** function o17excess(R17,R18) !use infotrac_phy, ONLY: use_iso USE isotopes_mod, ONLY: tnat,iso_O18,iso_O17 implicit none real R17,R18,o17excess if ((iso_O17.gt.0).and.(iso_O18.gt.0)) then o17excess=1e6*(log(R17/tnat(iso_o17)) & & -0.528*log(R18/tnat(iso_O18))) ! write(*,*) 'o17excess=',o17excess else write(*,*) 'iso_verif_egalite>deltaD 260: iso_O17.gt.0,18=', & & iso_O17.gt.0,iso_O18.gt.0 endif return end function o17excess ! **************** subroutine iso_verif_egalite_vect2D( & & xt,q,err_msg,ni,n,m) !use infotrac_phy, ONLY: use_iso USE isotopes_mod, ONLY: iso_eau implicit none ! inputs integer n,m,ni real q(n,m) real xt(ni,n,m) character*(*) err_msg ! locals integer iso_verif_egalite_nostop_loc integer i,j,ixt integer ifaux,jfaux !write(*,*) 'iso_verif_egalite_vect2D 1099 tmp: q(2,1),xt(iso_eau,2,1)=',q(2,1),xt(iso_eau,2,1) !write(*,*) 'ni,n,m=',ni,n,m,errmax,errmaxrel if (iso_eau.gt.0) then iso_verif_egalite_nostop_loc=0 do i=1,n do j=1,m if (abs(q(i,j)-xt(iso_eau,i,j)).gt.errmax) then if (abs((q(i,j)-xt(iso_eau,i,j))/ & & max(max(abs(q(i,j)),abs(xt(iso_eau,i,j))),1e-18)) & & .gt.errmaxrel) then iso_verif_egalite_nostop_loc=1 ifaux=i jfaux=j endif endif enddo !do j=1,m enddo !do i=1,n if (iso_verif_egalite_nostop_loc.eq.1) then write(*,*) 'erreur detectee par iso_verif_egalite_vect2D:' write(*,*) err_msg write(*,*) 'i,j=',ifaux,jfaux write(*,*) 'xt,q=',xt(iso_eau,ifaux,jfaux),q(ifaux,jfaux) stop endif endif #ifdef ISOVERIF call iso_verif_noNaN_vect2D(xt,err_msg,ni,n,m) #endif return end subroutine iso_verif_egalite_vect2D subroutine iso_verif_egalite_vect1D( & & xt,q,err_msg,ni,n) !use infotrac_phy, ONLY: use_iso USE isotopes_mod, ONLY: iso_eau implicit none ! inputs integer n,ni real q(n) real xt(ni,n) character*(*) err_msg ! locals integer iso_verif_egalite_nostop_loc integer i integer ifaux if (iso_eau.gt.0) then iso_verif_egalite_nostop_loc=0 do i=1,n if (abs(q(i)-xt(iso_eau,i)).gt.errmax) then if (abs((q(i)-xt(iso_eau,i))/ & & max(max(abs(q(i)),abs(xt(iso_eau,i))),1e-18)) & & .gt.errmaxrel) then iso_verif_egalite_nostop_loc=1 ifaux=i endif !if (abs((q(i)-xt(iso_eau,i))/ endif !if (abs(q(i)-xt(iso_eau,i)).gt.errmax) then enddo !do i=1,n if (iso_verif_egalite_nostop_loc.eq.1) then write(*,*) 'erreur detectee par iso_verif_egalite_vect2D:' write(*,*) err_msg write(*,*) 'i=',ifaux write(*,*) 'xt,q=',xt(iso_eau,ifaux),q(ifaux) stop endif !if (iso_verif_egalite_nostop.eq.1) then endif !if (iso_eau.gt.0) then end subroutine iso_verif_egalite_vect1D subroutine iso_verif_egalite_std_vect( & & a,b,err_msg,n,m,errmax,errmaxrel) implicit none ! inputs integer n,m,ni real a(n,m) real b(n,m) character*(*) err_msg real errmax,errmaxrel ! locals integer iso_verif_egalite_nostop_loc integer i,j integer ifaux,jfaux iso_verif_egalite_nostop_loc=0 do i=1,n do j=1,m if (abs(a(i,j)-b(i,j)).gt.errmax) then if (abs((a(i,j)-b(i,j))/ & & max(max(abs(a(i,j)),abs(b(i,j))),1e-18)) & & .gt.errmaxrel) then iso_verif_egalite_nostop_loc=1 ifaux=i jfaux=j endif endif enddo !do j=1,m enddo !do i=1,n if (iso_verif_egalite_nostop_loc.eq.1) then write(*,*) 'erreur detectee par iso_verif_egalite_vect2D:' write(*,*) err_msg write(*,*) 'i,j=',ifaux,jfaux write(*,*) 'a,b=',a(ifaux,jfaux),b(ifaux,jfaux) stop endif return end subroutine iso_verif_egalite_std_vect subroutine iso_verif_aberrant_vect2D( & & xt,q,err_msg,ni,n,m) !use infotrac_phy, ONLY: use_iso use isotopes_mod, ONLY: ridicule,tnat,iso_HDO implicit none ! inputs integer n,m,ni real q(n,m) real xt(ni,n,m) character*(*) err_msg ! locals integer iso_verif_aberrant_nostop_loc integer i,j integer ifaux,jfaux !real deltaD if (iso_HDO.gt.0) then iso_verif_aberrant_nostop_loc=0 do i=1,n do j=1,m if (q(i,j).gt.ridicule) then if (((xt(iso_HDO,i,j)/q(i,j)/tnat(iso_HDO)-1)*1000.0 & & .gt.deltalim).or. & & ((xt(iso_HDO,i,j)/q(i,j)/tnat(iso_HDO)-1)*1000.0 & & .lt.-borne)) then iso_verif_aberrant_nostop_loc=1 ifaux=i jfaux=j endif endif enddo !do j=1,m enddo !do i=1,n if (iso_verif_aberrant_nostop_loc.eq.1) then write(*,*) 'erreur detectee par iso_verif_aberrant_vect2D:' write(*,*) err_msg write(*,*) 'i,j=',ifaux,jfaux write(*,*) 'deltaD=',deltaD(xt(iso_HDO,ifaux,jfaux) & & /q(ifaux,jfaux)) write(*,*) 'xt(:,ifaux,jfaux)=',xt(:,ifaux,jfaux) stop endif endif !if (iso_HDO.gt.0) then end subroutine iso_verif_aberrant_vect2D subroutine iso_verif_aberrant_enc_vect2D( & & xt,q,err_msg,ni,n,m) !use infotrac_phy, ONLY: use_iso use isotopes_mod, ONLY: ridicule,tnat,iso_HDO implicit none ! inputs integer n,m,ni real q(n,m) real xt(ni,n,m) character*(*) err_msg ! locals integer iso_verif_aberrant_nostop_loc integer i,j integer ifaux,jfaux !real deltaD if (iso_HDO.gt.0) then iso_verif_aberrant_nostop_loc=0 do i=1,n do j=1,m if (q(i,j).gt.ridicule) then if (((xt(iso_HDO,i,j)/q(i,j)/tnat(iso_HDO)-1)*1000.0 & & .gt.deltalim).or. & & ((xt(iso_HDO,i,j)/q(i,j)/tnat(iso_HDO)-1)*1000.0 & & .lt.deltaDmin).or. & & (xt(iso_HDO,i,j).lt.-borne).or. & & (xt(iso_HDO,i,j).gt.borne)) then iso_verif_aberrant_nostop_loc=1 ifaux=i jfaux=j endif endif enddo !do j=1,m enddo !do i=1,n if (iso_verif_aberrant_nostop_loc.eq.1) then write(*,*) 'erreur detectee par ', & & 'iso_verif_aberrant_enc_vect2D:' write(*,*) err_msg write(*,*) 'i,j=',ifaux,jfaux write(*,*) 'deltaD=',deltaD(xt(iso_HDO,ifaux,jfaux) & & /q(ifaux,jfaux)) write(*,*) 'xt(:,ifaux,jfaux)=',xt(:,ifaux,jfaux) write(*,*) 'q(ifaux,jfaux)=',q(ifaux,jfaux) call abort_physic('isotopes_verif_mod','iso_verif_aberrant_enc_vect2D',1) endif endif !if (iso_HDO.gt.0) then end subroutine iso_verif_aberrant_enc_vect2D subroutine iso_verif_aberrant_enc_vect2D_ns( & & xt,q,err_msg,ni,n,m) !use infotrac_phy, ONLY: use_iso use isotopes_mod, ONLY: ridicule,tnat,iso_HDO implicit none ! inputs integer n,m,ni real q(n,m) real xt(ni,n,m) character*(*) err_msg ! locals integer iso_verif_aberrant_nostop_loc integer i,j integer ifaux,jfaux !real deltaD if (iso_HDO.gt.0) then iso_verif_aberrant_nostop_loc=0 do i=1,n do j=1,m if (q(i,j).gt.ridicule) then if (((xt(iso_HDO,i,j)/q(i,j)/tnat(iso_HDO)-1)*1000.0 & & .gt.deltalim).or. & & ((xt(iso_HDO,i,j)/q(i,j)/tnat(iso_HDO)-1)*1000.0 & & .lt.deltaDmin)) then iso_verif_aberrant_nostop_loc=1 ifaux=i jfaux=j endif endif enddo !do j=1,m enddo !do i=1,n if (iso_verif_aberrant_nostop_loc.eq.1) then write(*,*) 'erreur detectee par ', & & 'iso_verif_aberrant_vect2D_ns:' write(*,*) err_msg write(*,*) 'i,j=',ifaux,jfaux write(*,*) 'deltaD=',deltaD(xt(iso_HDO,ifaux,jfaux) & & /q(ifaux,jfaux)) write(*,*) 'xt(:,ifaux,jfaux)=',xt(:,ifaux,jfaux) ! stop endif endif !if (iso_HDO.gt.0) then end subroutine iso_verif_aberrant_enc_vect2D_ns subroutine iso_verif_aberrant_vect2Dch( & & xt,q,err_msg,ni,n,m,deltaDmax) !use infotrac_phy, ONLY: use_iso use isotopes_mod, ONLY: ridicule,tnat,iso_HDO implicit none ! inputs integer n,m,ni real q(n,m) real xt(ni,n,m) character*(*) err_msg real deltaDmax ! locals integer iso_verif_aberrant_nostop_loc integer i,j integer ifaux,jfaux !real deltaD if (iso_HDO.gt.0) then iso_verif_aberrant_nostop_loc=0 do i=1,n do j=1,m if (q(i,j).gt.ridicule) then if (((xt(iso_HDO,i,j)/q(i,j)/tnat(iso_HDO)-1)*1000.0 & & .gt.deltaDmax).or. & & ((xt(iso_HDO,i,j)/q(i,j)/tnat(iso_HDO)-1)*1000.0 & & .lt.-borne)) then iso_verif_aberrant_nostop_loc=1 ifaux=i jfaux=j endif endif enddo !do j=1,m enddo !do i=1,n if (iso_verif_aberrant_nostop_loc.eq.1) then write(*,*) 'erreur detectee par iso_verif_aberrant_vect2D:' write(*,*) err_msg write(*,*) 'i,j=',ifaux,jfaux write(*,*) 'deltaD=',deltaD(xt(iso_HDO,ifaux,jfaux) & & /q(ifaux,jfaux)) write(*,*) 'xt(:,ifaux,jfaux)=',xt(:,ifaux,jfaux) stop endif endif !if (iso_HDO.gt.0) then end subroutine iso_verif_aberrant_vect2Dch subroutine iso_verif_O18_aberrant_enc_vect2D( & & xt,q,err_msg,ni,n,m) !use infotrac_phy, ONLY: use_iso use isotopes_mod, ONLY: ridicule,tnat,iso_HDO,iso_O18 implicit none ! inputs integer n,m,ni real q(n,m) real xt(ni,n,m) character*(*) err_msg ! locals integer iso_verif_aberrant_nostop_loc integer i,j integer ifaux,jfaux real deltaDloc,deltaoloc,dexcessloc if ((iso_HDO.gt.0).and.(iso_O18.gt.0)) then iso_verif_aberrant_nostop_loc=0 do i=1,n do j=1,m if (q(i,j).gt.ridicule) then deltaDloc=(xt(iso_HDO,i,j)/q(i,j)/tnat(iso_hdo)-1)*1000 deltaoloc=(xt(iso_O18,i,j)/q(i,j)/tnat(iso_O18)-1)*1000 dexcessloc=deltaDloc-8*deltaoloc if ((deltaDloc.lt.deltaDmin).or.(deltaoloc.lt.deltaDmin/2.0).or. & & (deltaDloc.gt.deltalim).or.(deltaoloc.gt.deltalim/8.0).or. & & ((deltaDloc.gt.-500.0).and.((dexcessloc.lt.dexcess_min) & & .or.(dexcessloc.gt.dexcess_max)))) then iso_verif_aberrant_nostop_loc=1 ifaux=i jfaux=j write(*,*) 'deltaD,deltao,dexcess=',deltaDloc,deltaoloc,dexcessloc endif endif enddo !do j=1,m enddo !do i=1,n if (iso_verif_aberrant_nostop_loc.eq.1) then write(*,*) 'erreur detectee par ', & & 'iso_verif_aberrant_enc_vect2D:' write(*,*) err_msg write(*,*) 'i,j=',ifaux,jfaux write(*,*) 'xt(:,ifaux,jfaux)=',xt(:,ifaux,jfaux) write(*,*) 'q(ifaux,jfaux)=',q(ifaux,jfaux) call abort_physic('isotopes_verif_mod','iso_verif_aberrant_enc_vect2D',1) endif endif !if (iso_HDO.gt.0) then end subroutine iso_verif_O18_aberrant_enc_vect2D subroutine select_dim23_from4D(n1,n2,n3,n4, & & var,vec,i1,i4) implicit none ! inputs integer n1,n2,n3,n4 real var(n1,n2,n3,n4) integer i1,i4 ! outputs real vec(n2,n3) ! locals integer i2,i3 do i2=1,n2 do i3=1,n3 vec(i2,i3)=var(i1,i2,i3,i4) enddo enddo return end subroutine select_dim23_from4D subroutine select_dim4_from4D(ntime,nlev,nlat,nlon, & & var,vec,itime,ilev,ilat) implicit none ! inputs integer ntime,nlev,nlat,nlon real var(ntime,nlev,nlat,nlon) integer itime,ilev,ilat ! outputs real vec(nlon) ! locals integer ilon do ilon=1,nlon vec(ilon)=var(itime,ilev,ilat,ilon) enddo return end subroutine select_dim4_from4D subroutine select_dim5_from5D(n1,n2,n3,n4,n5, & & var,vec,i1,i2,i3,i4) implicit none ! inputs integer n1,n2,n3,n4,n5 real var(n1,n2,n3,n4,n5) integer i1,i2,i3,i4 ! outputs real vec(n5) ! locals integer i5 do i5=1,n5 vec(i5)=var(i1,i2,i3,i4,i5) enddo end subroutine select_dim5_from5D subroutine select_dim3_from3D(ntime,nlat,nlon, & & var,vec,itime,ilat) implicit none ! inputs integer ntime,nlat,nlon real var(ntime,nlat,nlon) integer itime,ilat ! outputs real vec(nlon) ! locals integer ilon do ilon=1,nlon vec(ilon)=var(itime,ilat,ilon) enddo end subroutine select_dim3_from3D subroutine select_dim23_from3D(n1,n2,n3, & & var,vec,i1) implicit none ! inputs integer n1,n2,n3 real var(n1,n2,n3) integer i1 ! outputs real vec(n2,n3) ! locals integer i2,i3 do i2=1,n2 do i3=1,n3 vec(i2,i3)=var(i1,i2,i3) enddo enddo end subroutine select_dim23_from3D subroutine putinto_dim23_from4D(n1,n2,n3,n4, & & var,vec,i1,i4) implicit none ! inputs integer n1,n2,n3,n4 real vec(n2,n3) integer i1,i4 ! inout real var(n1,n2,n3,n4) ! locals integer i2,i3 do i2=1,n2 do i3=1,n3 var(i1,i2,i3,i4)=vec(i2,i3) enddo enddo end subroutine putinto_dim23_from4D subroutine putinto_dim12_from4D(n1,n2,n3,n4, & & var,vec,i3,i4) implicit none ! inputs integer n1,n2,n3,n4 real vec(n1,n2) integer i3,i4 ! inout real var(n1,n2,n3,n4) ! locals integer i1,i2 do i1=1,n1 do i2=1,n2 var(i1,i2,i3,i4)=vec(i1,i2) enddo enddo end subroutine putinto_dim12_from4D subroutine putinto_dim23_from3D(n1,n2,n3, & & var,vec,i1) implicit none ! inputs integer n1,n2,n3 real vec(n2,n3) integer i1 ! inout real var(n1,n2,n3) ! locals integer i2,i3 do i2=1,n2 do i3=1,n3 var(i1,i2,i3)=vec(i2,i3) enddo enddo end subroutine putinto_dim23_from3D subroutine iso_verif_noNaN_par2D(x,err_msg,ni,n,m,ib,ie) implicit none ! si x est NaN, on affiche message ! d'erreur et return 1 si erreur ! input: integer n,m,ni,ib,ie real x(ni,n,m) character*(*) err_msg ! message d''erreur à afficher ! output ! locals integer i,j,ixt do i=ib,ie do j=1,m do ixt=1,ni if ((x(ixt,i,j).gt.-borne).and. & & (x(ixt,i,j).lt.borne)) then else !if ((x(ixt,i,j).gt.-borne).and. write(*,*) 'erreur detectee par iso_verif_nonNAN:' write(*,*) err_msg write(*,*) 'x,ixt,i,j=',x(ixt,i,j),ixt,i,j stop endif !if ((x(ixt,i,j).gt.-borne).and. enddo !do ixt=1,ni enddo !do j=1,m enddo !do i=1,n #ifdef ISOVERIF #else write(*,*) 'err_msg iso1772=',err_msg,': pourquoi verif?' stop #endif return end subroutine iso_verif_noNaN_par2D subroutine iso_verif_aberrant_enc_par2D( & & xt,q,err_msg,ni,n,m,ib,ie) !use infotrac_phy, ONLY: use_iso use isotopes_mod, ONLY: ridicule,tnat,iso_HDO implicit none ! inputs integer n,m,ni,ib,ie real q(n,m) real xt(ni,n,m) character*(*) err_msg ! locals integer iso_verif_aberrant_nostop_loc integer i,j integer ifaux,jfaux !real deltaD if (iso_HDO.gt.0) then iso_verif_aberrant_nostop_loc=0 do i=ib,ie do j=1,m if (q(i,j).gt.ridicule) then if (((xt(iso_HDO,i,j)/q(i,j)/tnat(iso_HDO)-1)*1000.0 & & .gt.deltalim).or. & & ((xt(iso_HDO,i,j)/q(i,j)/tnat(iso_HDO)-1)*1000.0 & & .lt.deltaDmin)) then iso_verif_aberrant_nostop_loc=1 ifaux=i jfaux=j endif endif enddo !do j=1,m enddo !do i=1,n if (iso_verif_aberrant_nostop_loc.eq.1) then write(*,*) 'erreur detectee par iso_verif_aberrant_par2D:' write(*,*) err_msg write(*,*) 'i,j=',ifaux,jfaux write(*,*) 'deltaD=',deltaD(xt(iso_HDO,ifaux,jfaux) & & /q(ifaux,jfaux)) write(*,*) 'xt(:,ifaux,jfaux)=',xt(:,ifaux,jfaux) write(*,*) 'q(ifaux,jfaux)=',q(ifaux,jfaux) stop endif endif !if (iso_HDO.gt.0) then end subroutine iso_verif_aberrant_enc_par2D subroutine iso_verif_egalite_par2D( & & xt,q,err_msg,ni,n,m,ib,ie) !use infotrac_phy, ONLY: use_iso USE isotopes_mod, ONLY: iso_eau implicit none ! inputs integer n,m,ni,ib,ie real q(n,m) real xt(ni,n,m) character*(*) err_msg ! locals integer iso_verif_egalite_nostop_loc integer i,j integer ifaux,jfaux if (iso_eau.gt.0) then iso_verif_egalite_nostop_loc=0 do i=ib,ie do j=1,m if (abs(q(i,j)-xt(iso_eau,i,j)).gt.errmax) then if (abs((q(i,j)-xt(iso_eau,i,j))/ & & max(max(abs(q(i,j)),abs(xt(iso_eau,i,j))),1e-18)) & & .gt.errmaxrel) then iso_verif_egalite_nostop_loc=1 ifaux=i jfaux=j endif endif enddo !do j=1,m enddo !do i=1,n if (iso_verif_egalite_nostop_loc.eq.1) then write(*,*) 'erreur detectee par iso_verif_egalite_vect2D:' write(*,*) err_msg write(*,*) 'i,j=',ifaux,jfaux write(*,*) 'xt,q=',xt(iso_eau,ifaux,jfaux),q(ifaux,jfaux) stop endif endif end subroutine iso_verif_egalite_par2D #ifdef ISOTRAC function iso_verif_traceur_choix_nostop(x,err_msg, & & errmax,errmaxrel,ridicule_trac,deltalimtrac) USE infotrac_phy, ONLY: ntraciso use isotopes_mod, ONLY: iso_HDO implicit none ! vérifier des choses sur les traceurs ! * toutes les zones donne t l'istope total ! * pas de deltaD aberrant ! inputs real x(ntraciso) character*(*) err_msg ! message d''erreur à afficher real errmax,errmaxrel,ridicule_trac,deltalimtrac ! output integer iso_verif_traceur_choix_nostop ! locals !integer iso_verif_traceur_noNaN_nostop !integer iso_verif_tracm_choix_nostop !integer iso_verif_tracdD_choix_nostop !integer iso_verif_tracpos_choix_nostop iso_verif_traceur_choix_nostop=0 ! verif noNaN if (iso_verif_traceur_noNaN_nostop(x,err_msg).eq.1) then iso_verif_traceur_choix_nostop=1 endif ! verif masse if (iso_verif_tracm_choix_nostop(x,err_msg, & & errmax,errmaxrel).eq.1) then iso_verif_traceur_choix_nostop=1 endif ! verif deltaD if (iso_HDO.gt.0) then if (iso_verif_tracdD_choix_nostop(x,err_msg, & & ridicule_trac,deltalimtrac).eq.1) then iso_verif_traceur_choix_nostop=1 endif endif !if (iso_HDO.gt.0) then ! verif pas aberramment negatif if (iso_verif_tracpos_choix_nostop(x,err_msg, & & 1e-5).eq.1) then iso_verif_traceur_choix_nostop=1 endif end function iso_verif_traceur_choix_nostop function iso_verif_tracnps_choix_nostop(x,err_msg, & & errmax,errmaxrel,ridicule_trac,deltalimtrac) USE infotrac_phy, ONLY: ntraciso USE isotopes_mod, ONLY: iso_HDO implicit none ! vérifier des choses sur les traceurs ! * toutes les zones donne t l'istope total ! * pas de deltaD aberrant ! on ne vérfie pas la positivité ! inputs real x(ntraciso) character*(*) err_msg ! message d''erreur à afficher real errmax,errmaxrel,ridicule_trac,deltalimtrac ! output integer iso_verif_tracnps_choix_nostop ! locals !integer iso_verif_traceur_noNaN_nostop !integer iso_verif_tracm_choix_nostop !integer iso_verif_tracdD_choix_nostop iso_verif_tracnps_choix_nostop=0 ! verif noNaN if (iso_verif_traceur_noNaN_nostop(x,err_msg).eq.1) then iso_verif_tracnps_choix_nostop=1 endif ! verif masse if (iso_verif_tracm_choix_nostop(x,err_msg, & & errmax,errmaxrel).eq.1) then iso_verif_tracnps_choix_nostop=1 endif ! verif deltaD if (iso_HDO.gt.0) then if (iso_verif_tracdD_choix_nostop(x,err_msg, & & ridicule_trac,deltalimtrac).eq.1) then iso_verif_tracnps_choix_nostop=1 endif endif ! if (iso_HDO.gt.0) then return end function iso_verif_tracnps_choix_nostop function iso_verif_tracpos_choix_nostop(x,err_msg,seuil) use infotrac_phy, ONLY: ntraciso,niso use isotrac_mod, only: index_iso,strtrac,index_zone use isotopes_mod, only: striso implicit none ! inputs real x(ntraciso) character*(*) err_msg ! message d''erreur à afficher real seuil ! output integer iso_verif_tracpos_choix_nostop ! locals integer lnblnk integer iiso,ixt !integer iso_verif_positif_choix_nostop iso_verif_tracpos_choix_nostop=0 do ixt=niso+1,ntraciso iiso=index_iso(ixt) if (iso_verif_positif_choix_nostop(x(ixt),seuil,err_msg// & & ', verif positif, iso'//striso(iiso) & & //strtrac(index_zone(ixt))).eq.1) then iso_verif_tracpos_choix_nostop=1 endif enddo end function iso_verif_tracpos_choix_nostop function iso_verif_traceur_noNaN_nostop(x,err_msg) use infotrac_phy, ONLY: ntraciso,niso use isotrac_mod, only: index_iso use isotopes_mod, only: striso implicit none ! on vérifie juste que pas NaN ! inputs real x(ntraciso) character*(*) err_msg ! message d''erreur à afficher ! output integer iso_verif_traceur_noNaN_nostop ! locals integer lnblnk integer iiso,ixt !integer iso_verif_nonaN_nostop iso_verif_traceur_noNaN_nostop=0 do ixt=niso+1,ntraciso iiso=index_iso(ixt) ! write(*,*) 'iso_verif_traceurs 154: iiso,ixt=',iiso,ixt if (iso_verif_noNaN_nostop(x(ixt),err_msg// & & ', verif trac no NaN, iso'//striso(iiso)) & & .eq.1) then iso_verif_traceur_noNaN_nostop=1 endif enddo end function iso_verif_traceur_noNaN_nostop function iso_verif_tracm_choix_nostop(x,err_msg, & & errmaxin,errmaxrelin) use infotrac_phy, ONLY: index_trac,ntraciso,niso use isotopes_mod, ONLY: ridicule,striso use isotrac_mod, only: ntraceurs_zone ! on vérifie juste bilan de masse implicit none ! inputs real x(ntraciso) character*(*) err_msg ! message d''erreur à afficher real errmaxin,errmaxrelin ! output integer iso_verif_tracm_choix_nostop ! locals !integer iso_verif_egalite_choix_nostop integer iiso,izone,ixt real xtractot iso_verif_tracm_choix_nostop=0 do iiso=1,niso xtractot=0.0 do izone=1,ntraceurs_zone ixt=index_trac(izone,iiso) xtractot=xtractot+x(ixt) enddo !do izone=1,ntraceurs_zone if (iso_verif_egalite_choix_nostop(xtractot,x(iiso), & & err_msg//', verif trac egalite, iso '//striso(iiso), & & errmaxin,errmaxrelin).eq.1) then write(*,*) 'iso_verif_traceur 202: x=',x ! write(*,*) 'xtractot=',xtractot iso_verif_tracm_choix_nostop=1 endif ! verif ajoutée le 19 fev 2011 if ((abs(xtractot).lt.ridicule**2).and. & & (abs(x(iiso)).gt.ridicule)) then write(*,*) err_msg,', verif masse traceurs, iso ', & & striso(iiso) write(*,*) 'iso_verif_traceur 209: x=',x ! iso_verif_tracm_choix_nostop=1 endif enddo !do iiso=1,ntraceurs_iso end function iso_verif_tracm_choix_nostop function iso_verif_tracdD_choix_nostop(x,err_msg, & & ridicule_trac,deltalimtrac) use infotrac_phy, ONLY: index_trac,ntraciso USE isotopes_mod, ONLY: iso_eau, iso_HDO use isotrac_mod, only: strtrac,ntraceurs_zone ! on vérifie juste deltaD implicit none ! inputs real x(ntraciso) character*(*) err_msg ! message d''erreur à afficher real ridicule_trac,deltalimtrac ! output integer iso_verif_tracdD_choix_nostop ! locals integer izone,ieau,ixt !integer iso_verif_aberrant_choix_nostop iso_verif_tracdD_choix_nostop=0 if ((iso_eau.gt.0).and.(iso_HDO.gt.0)) then do izone=1,ntraceurs_zone ieau=index_trac(izone,iso_eau) ixt=index_trac(izone,iso_HDO) if (iso_verif_aberrant_choix_nostop(x(ixt),x(ieau), & & ridicule_trac,deltalimtrac,err_msg// & & ', verif trac no aberrant zone '//strtrac(izone)) & & .eq.1) then iso_verif_tracdD_choix_nostop=1 endif ! if (x(ieau).gt.ridicule) then ! call iso_verif_aberrant(x(ixt)/x(ieau), ! : err_msg//', verif trac no aberrant zone ' ! : //strtrac(izone)) ! endif enddo !do izone=1,ntraceurs_zone endif ! if ((iso_eau.gt.0).and.(iso_HDO.gt.0)) then end function iso_verif_tracdD_choix_nostop subroutine iso_verif_trac17_q_deltaD(x,err_msg) use isotrac_mod, only: nzone_temp,option_traceurs USE infotrac_phy, ONLY: ntraciso implicit none ! inputs real x(ntraciso) character*(*) err_msg ! local integer iso_verif_tag17_q_deltaD_chns if ((option_traceurs.eq.17).or. & & (option_traceurs.eq.18)) then if (nzone_temp.ge.5) then if (iso_verif_tag17_q_deltaD_chns(x,err_msg).eq.1) then stop endif endif endif !if (option_traceurs.eq.17) then end subroutine iso_verif_trac17_q_deltaD subroutine iso_verif_traceur(x,err_msg) USE infotrac_phy, ONLY: ntraciso use isotrac_mod, only: ridicule_trac implicit none ! vérifier des choses sur les traceurs ! * toutes les zones donne t l'istope total ! * pas de deltaD aberrant ! on prend les valeurs pas défaut pour ! errmax,errmaxrel,ridicule_trac,deltalimtrac ! inputs real x(ntraciso) character*(*) err_msg ! message d''erreur à afficher ! locals !integer iso_verif_traceur_choix_nostop if (iso_verif_traceur_choix_nostop(x,err_msg, & & errmax,errmaxrel,ridicule_trac,deltalimtrac) & & .eq.1) then stop endif end subroutine iso_verif_traceur subroutine iso_verif_traceur_retourne3D(x,n1,n2,n3, & & i1,i2,i3,err_msg) USE infotrac_phy, ONLY: ntraciso use isotrac_mod, only: ridicule_trac implicit none ! vérifier des choses sur les traceurs ! * toutes les zones donne t l'istope total ! * pas de deltaD aberrant ! on prend les valeurs pas défaut pour ! errmax,errmaxrel,ridicule_trac,deltalimtrac ! inputs integer n1,n2,n3 real x(n1,n2,n3,ntraciso) character*(*) err_msg ! message d''erreur à afficher integer i1,i2,i3 ! locals !integer iso_verif_traceur_choix_nostop real xiso(ntraciso) call select_dim4_from4D(n1,n2,n3,ntraciso, & & x,xiso,i1,i2,i3) if (iso_verif_traceur_choix_nostop(xiso,err_msg, & & errmax,errmaxrel,ridicule_trac,deltalimtrac) & & .eq.1) then stop endif end subroutine iso_verif_traceur_retourne3D subroutine iso_verif_traceur_retourne4D(x,n1,n2,n3,n4, & & i1,i2,i3,i4,err_msg) USE infotrac_phy, ONLY: ntraciso use isotrac_mod, only: ridicule_trac implicit none ! vérifier des choses sur les traceurs ! * toutes les zones donne t l'istope total ! * pas de deltaD aberrant ! on prend les valeurs pas défaut pour ! errmax,errmaxrel,ridicule_trac,deltalimtrac ! inputs integer n1,n2,n3,n4 real x(n1,n2,n3,n4,ntraciso) character*(*) err_msg ! message d''erreur à afficher integer i1,i2,i3,i4 ! locals !integer iso_verif_traceur_choix_nostop real xiso(ntraciso) call select_dim5_from5D(n1,n2,n3,n4,ntraciso, & & x,xiso,i1,i2,i3,i4) if (iso_verif_traceur_choix_nostop(xiso,err_msg, & & errmax,errmaxrel,ridicule_trac,deltalimtrac) & & .eq.1) then stop endif end subroutine iso_verif_traceur_retourne4D subroutine iso_verif_traceur_retourne2D(x,n1,n2, & & i1,i2,err_msg) USE infotrac_phy, ONLY: ntraciso use isotrac_mod, only: ridicule_trac implicit none ! vérifier des choses sur les traceurs ! * toutes les zones donne t l'istope total ! * pas de deltaD aberrant ! on prend les valeurs pas défaut pour ! errmax,errmaxrel,ridicule_trac,deltalimtrac ! inputs integer n1,n2 real x(n1,n2,ntraciso) character*(*) err_msg ! message d''erreur à afficher integer i1,i2 ! locals !integer iso_verif_traceur_choix_nostop real xiso(ntraciso) call select_dim3_from3D(n1,n2,ntraciso, & & x,xiso,i1,i2) if (iso_verif_traceur_choix_nostop(xiso,err_msg, & & errmax,errmaxrel,ridicule_trac,deltalimtrac) & & .eq.1) then stop endif end subroutine iso_verif_traceur_retourne2D subroutine iso_verif_traceur_vect(x,n,m,err_msg) USE infotrac_phy, ONLY: ntraciso USE isotopes_mod, ONLY: iso_HDO implicit none ! vérifier des choses sur les traceurs ! * toutes les zones donne t l'istope total ! * pas de deltaD aberrant ! on prend les valeurs pas défaut pour ! errmax,errmaxrel,ridicule_trac,deltalimtrac ! inputs integer n,m real x(ntraciso,n,m) character*(*) err_msg ! message d''erreur à afficher ! locals logical iso_verif_traceur_nostop integer i,j,ixt,iiso,izone,ieau integer ifaux,jfaux,ixtfaux call iso_verif_traceur_noNaN_vect(x,n,m,err_msg) ! verif masse: iso_verif_tracm_choix_nostop call iso_verif_trac_masse_vect(x,n,m,err_msg,errmax,errmaxrel) ! verif deltaD: iso_verif_tracdD_choix_nostop if (iso_HDO.gt.0) then call iso_verif_tracdd_vect(x,n,m,err_msg) endif !if (iso_HDO.gt.0) then ! verif pas aberramment negatif: iso_verif_tracpos_choix_nostop call iso_verif_tracpos_vect(x,n,m,err_msg,1e-5) end subroutine iso_verif_traceur_vect subroutine iso_verif_tracnps_vect(x,n,m,err_msg) USE infotrac_phy, ONLY: ntraciso USE isotopes_mod, ONLY: iso_HDO implicit none ! vérifier des choses sur les traceurs ! * toutes les zones donne t l'istope total ! * pas de deltaD aberrant ! on prend les valeurs pas défaut pour ! errmax,errmaxrel,ridicule_trac,deltalimtrac ! inputs integer n,m real x(ntraciso,n,m) character*(*) err_msg ! message d''erreur à afficher ! locals logical iso_verif_traceur_nostop integer i,j,ixt,iiso,izone,ieau integer ifaux,jfaux,ixtfaux call iso_verif_traceur_noNaN_vect(x,n,m,err_msg) ! verif masse: iso_verif_tracm_choix_nostop call iso_verif_trac_masse_vect(x,n,m,err_msg,errmax,errmaxrel) ! verif deltaD: iso_verif_tracdD_choix_nostop if (iso_HDO.gt.0) then call iso_verif_tracdd_vect(x,n,m,err_msg) endif !if (iso_HDO.gt.0) then end subroutine iso_verif_tracnps_vect subroutine iso_verif_traceur_noNaN_vect(x,n,m,err_msg) USE infotrac_phy, ONLY: ntraciso,niso implicit none ! inputs integer n,m real x(ntraciso,n,m) character*(*) err_msg ! message d''erreur à afficher ! locals logical iso_verif_traceur_nostop integer i,j,ixt,iiso integer ifaux,jfaux,ixtfaux iso_verif_traceur_nostop=.false. ! verif noNaN: iso_verif_traceur_noNaN_nostop do j=1,m do i=1,n do ixt=niso+1,ntraciso if ((x(ixt,i,j).gt.-borne).and.(x(ixt,i,j).lt.borne)) then else !if ((x.gt.-borne).and.(x.lt.borne)) then iso_verif_traceur_nostop=.true. ifaux=i jfaux=i endif !if ((x.gt.-borne).and.(x.lt.borne)) then enddo !do ixt=niso+1,ntraciso enddo ! do i=1,n enddo ! do j=1,m if (iso_verif_traceur_nostop) then write(*,*) 'erreur detectée par iso_verif_nonNAN ', & & 'dans iso_verif_traceur_vect' write(*,*) '' write(*,*) err_msg write(*,*) 'x=',x(:,ifaux,jfaux) stop endif end subroutine iso_verif_traceur_noNaN_vect subroutine iso_verif_trac_masse_vect(x,n,m,err_msg, & & errmax,errmaxrel) USE infotrac_phy, ONLY: index_trac,ntraciso,niso use isotopes_mod, only: striso use isotrac_mod, only: ntraceurs_zone implicit none ! inputs integer n,m real x(ntraciso,n,m) character*(*) err_msg ! message d''erreur à afficher real errmax,errmaxrel ! locals logical iso_verif_traceur_nostop integer i,j,ixt,iiso,izone integer ifaux,jfaux,ixtfaux real xtractot(n,m) real xiiso(n,m) do iiso=1,niso do j=1,m do i=1,n xtractot(i,j)=0.0 xiiso(i,j)=x(iiso,i,j) do izone=1,ntraceurs_zone ixt=index_trac(izone,iiso) xtractot(i,j)=xtractot(i,j)+x(ixt,i,j) enddo !do izone=1,ntraceurs_zone enddo !do i=1,n enddo !do j=1,m call iso_verif_egalite_std_vect( & & xtractot,xiiso, & & err_msg//', verif trac egalite, iso '//striso(iiso), & & n,m,errmax,errmaxrel) enddo !do iiso=1,niso end subroutine iso_verif_trac_masse_vect subroutine iso_verif_tracdd_vect(x,n,m,err_msg) use infotrac_phy, only: index_trac,ntraciso,niso use isotopes_mod, only: iso_HDO,iso_eau use isotrac_mod, only: strtrac,ntraceurs_zone implicit none ! inputs integer n,m real x(ntraciso,n,m) character*(*) err_msg ! message d''erreur à afficher ! locals integer i,j,iiso,izone,ieau,ixt real xiiso(niso,n,m) real xeau(n,m) integer lnblnk if (iso_HDO.gt.0) then do izone=1,ntraceurs_zone ieau=index_trac(izone,iso_eau) do iiso=1,niso ixt=index_trac(izone,iiso) do j=1,m do i=1,n xiiso(iiso,i,j)=x(ixt,i,j) enddo !do i=1,n enddo !do j=1,m enddo !do iiso=1,niso do j=1,m do i=1,n xeau(i,j)=x(ieau,i,j) enddo !do i=1,n enddo !do j=1,m call iso_verif_aberrant_vect2Dch( & & xiiso,xeau,err_msg//strtrac(izone),niso,n,m, & & deltalimtrac) enddo !do izone=1,ntraceurs_zone endif !if (iso_HDO.gt.0) then end subroutine iso_verif_tracdd_vect subroutine iso_verif_tracpos_vect(x,n,m,err_msg,seuil) USE infotrac_phy, ONLY: ntraciso,niso implicit none ! inputs integer n,m real x(ntraciso,n,m) character*(*) err_msg ! message d''erreur à afficher real seuil ! locals integer i,j,ixt logical iso_verif_traceur_nostop integer ifaux,jfaux,ixtfaux iso_verif_traceur_nostop=.false. do j=1,m do i=1,n do ixt=niso+1,ntraciso if (x(ixt,i,j).lt.-seuil) then ifaux=i jfaux=j ixtfaux=ixt iso_verif_traceur_nostop=.true. endif enddo !do ixt=niso+1,ntraciso enddo !do i=1,n enddo !do j=1,m if (iso_verif_traceur_nostop) then write(*,*) 'erreur detectée par verif positif ', & & 'dans iso_verif_traceur_vect' write(*,*) '' write(*,*) err_msg write(*,*) 'x=',x(:,ifaux,jfaux) stop endif end subroutine iso_verif_tracpos_vect subroutine iso_verif_tracnps(x,err_msg) USE infotrac_phy, ONLY: ntraciso use isotrac_mod, only: ridicule_trac implicit none ! vérifier des choses sur les traceurs ! * toutes les zones donne t l'istope total ! * pas de deltaD aberrant ! on prend les valeurs pas défaut pour ! errmax,errmaxrel,ridicule_trac,deltalimtrac ! inputs real x(ntraciso) character*(*) err_msg ! message d''erreur à afficher ! locals !integer iso_verif_tracnps_choix_nostop if (iso_verif_tracnps_choix_nostop(x,err_msg, & & errmax,errmaxrel,ridicule_trac,deltalimtrac) & & .eq.1) then stop endif end subroutine iso_verif_tracnps subroutine iso_verif_tracpos_choix(x,err_msg,seuil) USE infotrac_phy, ONLY: ntraciso implicit none ! vérifier des choses sur les traceurs ! * toutes les zones donne t l'istope total ! * pas de deltaD aberrant ! on prend les valeurs pas défaut pour ! errmax,errmaxrel,ridicule_trac,deltalimtrac ! inputs real x(ntraciso) character*(*) err_msg ! message d''erreur à afficher real seuil ! locals !integer iso_verif_tracpos_choix_nostop if (iso_verif_tracpos_choix_nostop(x,err_msg,seuil) & & .eq.1) then stop endif end subroutine iso_verif_tracpos_choix subroutine iso_verif_traceur_choix(x,err_msg, & & errmax,errmaxrel,ridicule_trac_loc,deltalimtrac) USE infotrac_phy, ONLY: ntraciso implicit none ! vérifier des choses sur les traceurs ! * toutes les zones donne t l'istope total ! * pas de deltaD aberrant ! inputs real x(ntraciso) character*(*) err_msg ! message d''erreur à afficher real errmax,errmaxrel,ridicule_trac_loc,deltalimtrac ! locals !integer iso_verif_traceur_choix_nostop if (iso_verif_traceur_choix_nostop(x,err_msg, & & errmax,errmaxrel,ridicule_trac_loc,deltalimtrac) & & .eq.1) then stop endif end subroutine iso_verif_traceur_choix function iso_verif_traceur_nostop(x,err_msg) USE infotrac_phy, ONLY: ntraciso use isotrac_mod, only: ridicule_trac !use isotopes_verif, only: errmax,errmaxrel,deltalimtrac implicit none ! vérifier des choses sur les traceurs ! * toutes les zones donne t l'istope total ! * pas de deltaD aberrant ! on prend les valeurs pas défaut pour ! errmax,errmaxrel,ridicule_trac,deltalimtrac ! inputs real x(ntraciso) character*(*) err_msg ! message d''erreur à afficher ! output integer iso_verif_traceur_nostop ! locals !integer iso_verif_traceur_choix_nostop iso_verif_traceur_nostop= & & iso_verif_traceur_choix_nostop(x,err_msg, & & errmax,errmaxrel,ridicule_trac,deltalimtrac) end function iso_verif_traceur_nostop subroutine iso_verif_traceur_justmass(x,err_msg) USE infotrac_phy, ONLY: ntraciso implicit none ! on vérifie que noNaN et masse ! on prend les valeurs pas défaut pour ! errmax,errmaxrel,ridicule_trac,deltalimtrac ! inputs real x(ntraciso) character*(*) err_msg ! message d''erreur à afficher ! locals !integer iso_verif_traceur_noNaN_nostop !integer iso_verif_tracm_choix_nostop ! verif noNaN if (iso_verif_traceur_noNaN_nostop(x,err_msg).eq.1) then stop endif ! verif masse if (iso_verif_tracm_choix_nostop(x,err_msg, & & errmax,errmaxrel).eq.1) then stop endif end subroutine iso_verif_traceur_justmass function iso_verif_traceur_jm_nostop(x,err_msg) USE infotrac_phy, ONLY: ntraciso implicit none ! on vérifie que noNaN et masse ! on prend les valeurs pas défaut pour ! errmax,errmaxrel,ridicule_trac,deltalimtrac ! inputs real x(ntraciso) character*(*) err_msg ! message d''erreur à afficher ! output integer iso_verif_traceur_jm_nostop ! locals ! integer iso_verif_traceur_noNaN_nostop !integer iso_verif_tracm_choix_nostop iso_verif_traceur_jm_nostop=0 ! ! verif noNaN ! if (iso_verif_traceur_noNaN_nostop(x,err_msg).eq.1) then ! iso_verif_traceur_jm_nostop=1 ! endif ! verif masse if (iso_verif_tracm_choix_nostop(x,err_msg, & & errmax,errmaxrel).eq.1) then iso_verif_traceur_jm_nostop=1 endif end function iso_verif_traceur_jm_nostop function iso_verif_tag17_q_deltaD_chns(x,err_msg) USE infotrac_phy, ONLY: index_trac,ntraciso use isotopes_mod, ONLY: iso_HDO,iso_eau,ridicule use isotrac_mod, only: nzone_temp,option_traceurs implicit none ! inputs real x(ntraciso) character*(*) err_msg ! output integer iso_verif_tag17_q_deltaD_chns ! locals !integer iso_verif_positif_nostop !real deltaD integer ieau,ixt,ieau1 iso_verif_tag17_q_deltaD_chns=0 if ((option_traceurs.eq.17).or. & & (option_traceurs.eq.18)) then ! verifier que deltaD du tag de la couche la plus haute < ! 200 permil, et vérifier que son q est inférieur à ieau=index_trac(nzone_temp,iso_eau) ixt=index_trac(nzone_temp,iso_HDO) if (x(ieau).gt.ridicule) then if (iso_verif_positif_nostop(-200.0-deltaD(x(ixt)/x(ieau)), & & err_msg//': deltaDt05 trop fort').eq.1) then write(*,*) 'x=',x iso_verif_tag17_q_deltaD_chns=1 endif !if (iso_verif_positif_nostop(-200.0-deltaD(x(ixt),x(ieau)), endif !if (x(ieau).gt.ridicule) then if (iso_verif_positif_nostop(2.0e-3-x(ieau), & & err_msg//': qt05 trop fort').eq.1) then write(*,*) 'x=',x iso_verif_tag17_q_deltaD_chns=1 endif !if (iso_verif_positif_nostop(1.0e-3-x(ieau), ! on vérifie que si q est petit, alors qt01 fait moins de 10% if (x(iso_eau).lt.2.0e-3) then ieau1= index_trac(1,iso_eau) if (iso_verif_positif_nostop(0.1-(x(ieau1)/x(iso_eau)), & & err_msg//': qt01 trop abondant').eq.1) then write(*,*) 'x=',x iso_verif_tag17_q_deltaD_chns=1 endif ! if (iso_verif_positif(0.1-(x(ixt)/x(ieau)), endif !if (x(ieau).lt.2.0e-3) then endif !if (option_traceurs.eq.17) then end function iso_verif_tag17_q_deltaD_chns subroutine iso_verif_tag17_q_deltaD_vect(x,n,m,err_msg) USE infotrac_phy, ONLY: index_trac,ntraciso USE isotopes_mod, ONLY: tnat,iso_eau, ridicule,iso_HDO use isotrac_mod, only: option_traceurs,nzone_temp implicit none ! inputs integer n,m real x(ntraciso,n,m) character*(*) err_msg ! locals !integer iso_verif_positif_nostop !real deltaD integer ieau,ixt,ieau1 integer i,k if ((option_traceurs.eq.17).or. & & (option_traceurs.eq.18)) then ! verifier que deltaD du tag de la couche la plus haute < ! 200 permil, et vérifier que son q est inférieur à ieau=index_trac(nzone_temp,iso_eau) ixt=index_trac(nzone_temp,iso_HDO) ieau1=index_trac(1,iso_eau) do i=1,n do k=1,m if (x(ieau,i,k).gt.ridicule) then if ((x(ixt,i,k)/x(ieau,i,k)/tnat(iso_HDO)-1)*1000 & & .gt.-200.0) then write(*,*) err_msg,', vect:deltaDt05 trop fort' write(*,*) 'i,k=',i,k write(*,*) 'x(:,i,k)=',x(:,i,k) stop endif !if (iso_verif_positif_nostop(-200.0-deltaD(x(ixt),x(ieau)), endif !if (x(ieau).gt.ridicule) then if (x(ieau,i,k).gt.2.0e-3) then write(*,*) err_msg,', vect:qt05 trop fort' write(*,*) 'i,k=',i,k write(*,*) 'x(:,i,k)=',x(:,i,k) stop endif !if (iso_verif_positif_nostop(1.0e-3-x(ieau), if (x(iso_eau,i,k).lt.2.0e-3) then if (x(ieau1,i,k)/x(iso_eau,i,k).gt.0.1) then write(*,*) err_msg,', vect: qt01 trop abondant' write(*,*) 'i,k=',i,k write(*,*) 'ieau1,iso_eau,x(ieau1,i,k),', & & 'x(iso_eau,i,k)=',ieau1,iso_eau, & & x(ieau1,i,k),x(iso_eau,i,k) write(*,*) 'x(:,i,k)=',x(:,i,k) stop endif !if (x(ieau1,i,k)/x(iso_eau,i,k).gt.0.1) then endif enddo !do k=1,m enddo !do i=1,n endif !if (option_traceurs.eq.17) then end subroutine iso_verif_tag17_q_deltaD_vect subroutine iso_verif_tag17_q_deltaD_vect_ret3D(x,n,m,nq,err_msg) USE infotrac_phy, ONLY: index_trac,ntraciso USE isotopes_mod, ONLY: tnat,iso_eau,iso_HDO,ridicule use isotrac_mod, only: option_traceurs,nzone_temp implicit none ! inputs integer n,m,nq real x(n,m,nq,ntraciso) character*(*) err_msg ! locals !integer iso_verif_positif_nostop !real deltaD integer ieau,ixt,ieau1 integer i,k,iq if ((option_traceurs.eq.17).or. & & (option_traceurs.eq.18)) then ! verifier que deltaD du tag de la couche la plus haute < ! 200 permil, et vérifier que son q est inférieur à ieau=index_trac(nzone_temp,iso_eau) ixt=index_trac(nzone_temp,iso_HDO) ieau1=index_trac(1,iso_eau) do iq=1,nq do i=1,n do k=1,m if (x(i,k,iq,ieau).gt.ridicule) then if ((x(i,k,iq,ixt)/x(i,k,iq,ieau)/tnat(iso_HDO)-1)*1000 & & .gt.-200.0) then write(*,*) err_msg,', vect:deltaDt05 trop fort' write(*,*) 'i,k=',i,k write(*,*) 'x(i,k,iq,:)=',x(i,k,iq,:) stop endif !if (iso_verif_positif_nostop(-200.0-deltaD(x(ixt),x(ieau)), endif !if (x(ieau).gt.ridicule) then if (x(i,k,iq,ieau).gt.2.0e-3) then write(*,*) err_msg,', vect:qt05 trop fort' write(*,*) 'i,k=',i,k write(*,*) 'x(i,k,iq,:)=',x(i,k,iq,:) stop endif !if (iso_verif_positif_nostop(1.0e-3-x(ieau), if (x(i,k,iq,iso_eau).lt.2.0e-3) then if (x(i,k,iq,ieau1)/x(i,k,iq,iso_eau).gt.0.1) then write(*,*) err_msg,', vect: qt01 trop abondant' write(*,*) 'i,k=',i,k write(*,*) 'ieau1,iso_eau,x(i,k,iq,ieau1),', & & 'x(i,k,iq,ieau)=',ieau1,iso_eau, & & x(i,k,iq,ieau1),x(i,k,iq,iso_eau) write(*,*) 'x(i,k,iq,:)=',x(i,k,iq,:) stop endif !if (x(ieau1,i,k)/x(iso_eau,i,k).gt.0.1) then endif enddo !do k=1,m enddo !do i=1,n enddo ! do iq=1,nq endif !if (option_traceurs.eq.17) then end subroutine iso_verif_tag17_q_deltaD_vect_ret3D #endif ! endif ISOTRAC END MODULE isotopes_verif_mod #endif ! endif ISOVERIF