ANALYSE
NAME
ANALYSE
DESCRIPTION
Package data : This package implements all analysis availabale on QTLMap. +--------+ | | +--------+
m_qtlmap_analyse
[ Top ] [ ANALYSE ] [ Modules ]
NAME
m_qtlmap_analyse -- Interface analysis module using by the main program.
SYNOPSIS
DESCRIPTION
NOTES
SEE ALSO
analyse
[ Top ] [ m_qtlmap_analyse ] [ Subroutines ]
NAME
analyse
DESCRIPTION
Main function to analyse a user dataset with a specific analysis.
INPUTS
opt_calcul : analysis id opt_qtl : number of qtl hypothesis (analysis can be not compatible with this option) hdam : LD option, take dam haplotypes in considration opt_sim : contexte of the analysis ( COMMON_ANALYSE | SIMULATION | TRANSCIPTOME_ANALYSE )
OUTPUTS
lrtsol : likelihood ratio test information (see TYPE_LRT_SOLUTION) listincsol : solution for each hypothesis finded by the analysis (see TYPE_INCIDENCE_SOLUTION) rhoi : residual correlation matrix finded for each hypothesis
SOURCE
237 ! subroutine analyse(opt_calcul,opt_qtl,lrtsol,listincsol,rhoi,hdam,opt_sim,nopoly) 238 subroutine analyse(opt_calcul,opt_qtl,lrtsol,listincsol,rhoi,hdam,opt_sim) 239 use m_qtlmap_incidence 240 integer , intent(in) :: opt_calcul,opt_qtl 241 integer ,intent(in) :: opt_sim 242 type(TYPE_LRT_SOLUTION) , intent(inout) ,dimension(ncar,opt_qtl) :: lrtsol 243 type(TYPE_INCIDENCE_SOLUTION) ,intent(inout) ,dimension(ncar,opt_qtl+1) :: listincsol 244 real (kind=dp) ,dimension(size(carac),size(carac)), intent(out) :: rhoi 245 ! logical ,intent(in) :: hdam,nopoly 246 logical ,intent(in) :: hdam 247 integer :: ic 248 logical :: car_real 249 250 rhoi=0.d0 251 252 if ( opt_qtl <= 0) then 253 call stop_application("bad definition of opt_qtl ["//trim(str(opt_qtl))//"]") 254 end if 255 256 ! arret si l'utilisateur demande d'analyser que des variables enti�res ou cat�gorielles selon 257 ! une option non disponible (opt_qtl=2, opt_calcul=ANALYSE_MULTITRAIT ou ANALYSE_UNITRAIT_MODLIN_cox 258 259 car_real=.false. 260 do ic=1,size(carac) 261 if((natureY(ic) == 'r') .or.( natureY(ic) == 'a')) car_real=.true. 262 end do 263 if(.not.car_real) then 264 if (opt_calcul == ANALYSE_UNITRAIT_MODLIN_cox .or.& 265 opt_calcul == ANALYSE_MULTITRAIT) then 266 call stop_application("Analysis option not available [OPT_CALCUL="//trim(str(opt_calcul))//"]"& 267 //"[OPT_QTL="//trim(str(opt_qtl))//"] ") 268 end if 269 end if 270 271 if ( opt_calcul == ANALYSE_UNITRAIT .and. opt_qtl == 1) then 272 call opti_unitrait(lrtsol,listincsol,opt_sim,1) 273 else if ( opt_calcul == ANALYSE_UNITRAIT .and. opt_qtl == 2 ) then 274 call opti_unitrait(lrtsol,listincsol,opt_sim,2) 275 else if ( opt_calcul == ANALYSE_UNITRAIT_MODLIN .and. opt_qtl == 1) then 276 call opti_unitrait_modlin(lrtsol,listincsol,opt_sim) 277 if (opt_sim /= TRANSCIPTOME_ANALYSE) call opti_unitrait_discret(lrtsol,listincsol,opt_sim) 278 else if ( opt_calcul == ANALYSE_UNITRAIT_LDLA .and. opt_qtl == 1) then 279 call opti_unitrait_modlin_new(lrtsol(:,1),listincsol,"hete","LDLA",hdam,opt_sim) 280 else if ( opt_calcul == ANALYSE_UNITRAIT_LD .and. opt_qtl == 1) then 281 call opti_unitrait_modlin_new(lrtsol(:,1),listincsol,"hete","LD ",hdam,opt_sim) 282 else if ( opt_calcul == ANALYSE_UNITRAIT_LA .and. opt_qtl == 1) then 283 call opti_unitrait_modlin_new(lrtsol(:,1),listincsol,"hete","LA ",hdam,opt_sim) 284 else if ( opt_calcul == ANALYSE_UNITRAIT_LDJH .and. opt_qtl == 1) then 285 call opti_unitrait_modlin_new(lrtsol(:,1),listincsol,"hete","LDJH",hdam,opt_sim) 286 else if ( opt_calcul == ANALYSE_UNITRAIT_LINEAR_LA_HOMO .and. opt_qtl >= 1) then 287 call opti_unitrait_incidence& 288 ! (lrtsol,listincsol,opt_sim,INCIDENCE_TYPE_CLASSIC,LINEAR_HOMOSCEDASTIC,opt_qtl,hdam,nopoly) 289 (1,ncar,lrtsol,listincsol,opt_sim,INCIDENCE_TYPE_CLASSIC,LINEAR_HOMOSCEDASTIC,opt_qtl,hdam) 290 else if ( opt_calcul == ANALYSE_UNITRAIT_LINEAR_LA_HETERO .and. opt_qtl >= 1) then 291 call opti_unitrait_incidence& 292 ! (lrtsol,listincsol,opt_sim,INCIDENCE_TYPE_CLASSIC,LINEAR_HETEROSCEDASTIC,opt_qtl,hdam,nopoly) 293 (1,ncar,lrtsol,listincsol,opt_sim,INCIDENCE_TYPE_CLASSIC,LINEAR_HETEROSCEDASTIC,opt_qtl,hdam) 294 else if ( opt_calcul == ANALYSE_UNITRAIT_LINEAR_LD_HOMO .and. opt_qtl >= 1) then 295 call opti_unitrait_incidence& 296 ! (lrtsol,listincsol,opt_sim,INCIDENCE_TYPE_LD,LINEAR_HOMOSCEDASTIC,opt_qtl,hdam,nopoly) 297 (1,ncar,lrtsol,listincsol,opt_sim,INCIDENCE_TYPE_LD,LINEAR_HOMOSCEDASTIC,opt_qtl,hdam) 298 else if ( opt_calcul == ANALYSE_UNITRAIT_LINEAR_LD_HETERO .and. opt_qtl >= 1) then 299 call opti_unitrait_incidence& 300 ! (lrtsol,listincsol,opt_sim,INCIDENCE_TYPE_LD,LINEAR_HETEROSCEDASTIC,opt_qtl,hdam,nopoly) 301 (1,ncar,lrtsol,listincsol,opt_sim,INCIDENCE_TYPE_LD,LINEAR_HETEROSCEDASTIC,opt_qtl,hdam) 302 else if ( opt_calcul == ANALYSE_UNITRAIT_LINEAR_LDLA_HOMO .and. opt_qtl >= 1) then 303 call opti_unitrait_incidence& 304 ! (lrtsol,listincsol,opt_sim,INCIDENCE_TYPE_LDLA,LINEAR_HOMOSCEDASTIC,opt_qtl,hdam,nopoly) 305 (1,ncar,lrtsol,listincsol,opt_sim,INCIDENCE_TYPE_LDLA,LINEAR_HOMOSCEDASTIC,opt_qtl,hdam) 306 else if ( opt_calcul == ANALYSE_UNITRAIT_LINEAR_LDLA_HETERO .and. opt_qtl >= 1) then 307 call opti_unitrait_incidence& 308 ! (lrtsol,listincsol,opt_sim,INCIDENCE_TYPE_LDLA,LINEAR_HETEROSCEDASTIC,opt_qtl,hdam,nopoly) 309 (1,ncar,lrtsol,listincsol,opt_sim,INCIDENCE_TYPE_LDLA,LINEAR_HETEROSCEDASTIC,opt_qtl,hdam) 310 else if ( opt_calcul == ANALYSE_UNITRAIT_MODLIN_cox .and. opt_qtl == 1) then 311 call opti_unitrait_modlin_cox(lrtsol,listincsol,opt_sim) 312 else if ( opt_calcul == ANALYSE_MULTITRAIT .and. opt_qtl == 1) then 313 call opti_multitraits(lrtsol(1,1),listincsol,rhoi,opt_sim) 314 else if ( opt_calcul == ANALYSE_MULTITRAIT_DA .and. opt_qtl == 1) then 315 call opti_multitrait_DA(lrtsol(1,1),opt_sim) 316 else if ( opt_calcul == ANALYSE_2QTL_INTERACTION .and. opt_qtl == 2) then 317 call opti_unitrait_interaction(lrtsol,opt_sim) 318 else if ( opt_calcul == ANALYSE_UNITRAIT_CONTINGENCE .and. opt_qtl >= 1) then 319 call opti_unitrait_incidence& 320 ! (lrtsol,listincsol,opt_sim,INCIDENCE_TYPE_CLASSIC,OPTIM_HETEROSCEDASTIC,opt_qtl,hdam,nopoly) 321 (1,ncar,lrtsol,listincsol,opt_sim,INCIDENCE_TYPE_CLASSIC,OPTIM_HETEROSCEDASTIC,opt_qtl,hdam) 322 else if ( opt_calcul == ANALYSE_TRAIT_COV_CONTINGENCE .and. opt_qtl >= 1) then 323 call opti_unitrait_incidence& 324 ! (lrtsol,listincsol,opt_sim,INCIDENCE_TRAIT_COV,OPTIM_HETEROSCEDASTIC,opt_qtl,hdam,nopoly) 325 (1,ncar,lrtsol,listincsol,opt_sim,INCIDENCE_TRAIT_COV,OPTIM_HETEROSCEDASTIC,opt_qtl,hdam) 326 else if ( opt_calcul == ANALYSE_MULTITRAIT_INCIDENCE .and. opt_qtl >= 1) then 327 call opti_multitrait_modlin(lrtsol(1,:),listincsol,opt_sim,opt_qtl,ANALYSE_MULTITRAIT_INCIDENCE) 328 else if ( opt_calcul == ANALYSE_MULTITRAIT_INCIDENCE_LU .and. opt_qtl >= 1) then 329 call opti_multitrait_modlin(lrtsol(1,:),listincsol,opt_sim,opt_qtl,ANALYSE_MULTITRAIT_INCIDENCE_LU) 330 else if ( opt_calcul == ANALYSE_DEV_1 .and. opt_qtl == 1 ) then 331 call opti_unitrait_dev1(lrtsol,listincsol,opt_sim,OPTIM_HETEROSCEDASTIC,opt_qtl,hdam) 332 else 333 call stop_application("Analysis option not available [OPT_CALCUL="//trim(str(opt_calcul))//"]"& 334 //"[OPT_QTL="//trim(str(opt_qtl))//"]") 335 end if 336 337 if (opt_sim /= SIMULATION) call print_allelic_origin 338 339 340 end subroutine analyse
analyse_cuda
[ Top ] [ m_qtlmap_analyse ] [ Subroutines ]
NAME
analyse_cuda
DESCRIPTION
Main function to analyse a user dataset with a specific analysis.
INPUTS
opt_calcul : analysis id opt_qtl : number of qtl hypothesis (analysis can be not compatible with this option) hdam : LD option, take dam haplotypes in considration opt_sim : contexte of the analysis ( COMMON_ANALYSE | SIMULATION | TRANSCIPTOME_ANALYSE )
OUTPUTS
lrtsol : likelihood ratio test information (see TYPE_LRT_SOLUTION) listincsol : solution for each hypothesis finded by the analysis (see TYPE_INCIDENCE_SOLUTION) rhoi : residual correlation matrix finded for each hypothesis
SOURCE
360 #ifdef MANAGE_CUDA 361 subroutine analyse_cuda(nsim,opt_sim,opt_calcul,opt_qtl,lrtsol,listincsol,listLrtSolSim,ySIMUL,hdam) 362 use m_qtlmap_incidence 363 integer , intent(in) :: nsim,opt_calcul,opt_sim,opt_qtl 364 type(TYPE_LRT_SOLUTION) , intent(inout) ,dimension(ncar,opt_qtl) :: lrtsol 365 type(TYPE_INCIDENCE_SOLUTION) ,intent(inout) ,dimension(ncar,opt_qtl+1) :: listincsol 366 type(TYPE_LRT_SOLUTION) , intent(inout) ,dimension(nsim,ncar,opt_qtl) :: listLrtSolSim 367 real (kind=dp), dimension (:,:,:),allocatable , intent(in) :: ySIMUL 368 logical ,intent(in) :: hdam 369 integer :: ic 370 371 if ( opt_qtl <= 0) then 372 call stop_application("bad definition of opt_qtl ["//trim(str(opt_qtl))//"]") 373 end if 374 375 376 if ( opt_calcul == ANALYSE_UNITRAIT_LINEAR_LA_HOMO .and. opt_qtl >= 1) then 377 call opti_unitrait_incidence_cuda& 378 (nsim,lrtsol,listincsol,listLrtSolSim,opt_sim,INCIDENCE_TYPE_CLASSIC,LINEAR_HOMOSCEDASTIC,opt_qtl,ySIMUL,hdam) 379 else if ( opt_calcul == ANALYSE_UNITRAIT_LINEAR_LA_HETERO .and. opt_qtl >= 1) then 380 call opti_unitrait_incidence_cuda& 381 (nsim,lrtsol,listincsol,listLrtSolSim,opt_sim,INCIDENCE_TYPE_CLASSIC,LINEAR_HETEROSCEDASTIC,opt_qtl,ySIMUL,hdam) 382 else if ( opt_calcul == ANALYSE_UNITRAIT_LINEAR_LD_HOMO .and. opt_qtl >= 1) then 383 call opti_unitrait_incidence_cuda& 384 (nsim,lrtsol,listincsol,listLrtSolSim,opt_sim,INCIDENCE_TYPE_LD,LINEAR_HOMOSCEDASTIC,opt_qtl,ySIMUL,hdam) 385 else if ( opt_calcul == ANALYSE_UNITRAIT_LINEAR_LD_HETERO .and. opt_qtl >= 1) then 386 call opti_unitrait_incidence_cuda& 387 (nsim,lrtsol,listincsol,listLrtSolSim,opt_sim,INCIDENCE_TYPE_LD,LINEAR_HETEROSCEDASTIC,opt_qtl,ySIMUL,hdam) 388 else if ( opt_calcul == ANALYSE_UNITRAIT_LINEAR_LDLA_HOMO .and. opt_qtl >= 1) then 389 call opti_unitrait_incidence_cuda& 390 (nsim,lrtsol,listincsol,listLrtSolSim,opt_sim,INCIDENCE_TYPE_LDLA,LINEAR_HOMOSCEDASTIC,opt_qtl,ySIMUL,hdam) 391 else if ( opt_calcul == ANALYSE_UNITRAIT_LINEAR_LDLA_HETERO .and. opt_qtl >= 1) then 392 call opti_unitrait_incidence_cuda& 393 (nsim,lrtsol,listincsol,listLrtSolSim,opt_sim,INCIDENCE_TYPE_LDLA,LINEAR_HETEROSCEDASTIC,opt_qtl,ySIMUL,hdam) 394 else 395 call stop_application("Analysis option not available [OPT_CALCUL="//trim(str(opt_calcul))//"]"& 396 //"[OPT_QTL="//trim(str(opt_qtl))//"]") 397 end if 398 399 call print_allelic_origin 400 401 402 end subroutine analyse_cuda 403 #endif
get_name_analyse
[ Top ] [ m_qtlmap_analyse ] [ Subroutines ]
NAME
get_name_analyse
DESCRIPTION
Returns the string name of a analysis identified by the calcul id.
INPUTS
opt_calcul : analysis id
RETURN
name of the analysis
SOURCE
70 function get_name_analyse(opt_calcul) result (name) 71 integer , intent(in) :: opt_calcul 72 character(len=300) :: name 73 SELECT CASE (opt_calcul) 74 CASE (ANALYSE_UNITRAIT) 75 name ="UNITRAIT ANALYSIS" 76 CASE (ANALYSE_UNITRAIT_MODLIN) 77 name = "MODLIN ANALYSIS" 78 CASE (ANALYSE_MULTITRAIT) 79 name = "MULTI ANALYSIS" 80 CASE (ANALYSE_UNITRAIT_LINEAR_LA_HOMO) 81 name = "LINEAR LA HOMOSCEDASTIC" 82 CASE (ANALYSE_UNITRAIT_LINEAR_LA_HETERO) 83 name = "LINEAR LA HETEROSCEDASTIC" 84 CASE (ANALYSE_UNITRAIT_LINEAR_LD_HOMO) 85 name = "LINEAR LD HOMOSCEDASTIC" 86 CASE (ANALYSE_UNITRAIT_LINEAR_LD_HETERO) 87 name = "LINEAR LD HETEROSCEDASTIC" 88 CASE (ANALYSE_UNITRAIT_LINEAR_LDLA_HOMO) 89 name = "LINEAR LDLA HOMOSCEDASTIC" 90 CASE (ANALYSE_UNITRAIT_LINEAR_LDLA_HETERO) 91 name = "LINEAR LDLA HETEROSCEDASTIC" 92 CASE( ANALYSE_UNITRAIT_LDLA) 93 name = "LDLA ANALYSIS" 94 CASE( ANALYSE_UNITRAIT_LDJH) 95 name = "LDJH ANALYSIS" 96 CASE( ANALYSE_UNITRAIT_LD) 97 name = "LD ANALYSIS" 98 CASE( ANALYSE_UNITRAIT_LA) 99 name = "LA ANALYSIS (MODLIN)" 100 CASE (ANALYSE_UNITRAIT_MODLIN_COX) 101 name = "COX ANALYSIS" 102 CASE (ANALYSE_MULTITRAIT_DA) 103 name = "DISCRIMINANT ANALYSIS" 104 CASE (ANALYSE_2QTL_INTERACTION) 105 name = "2QTL INTERACTION ANALYSIS" 106 CASE (ANALYSE_UNITRAIT_CONTINGENCE) 107 name = "NEW MODLIN" 108 CASE (ANALYSE_TRAIT_COV_CONTINGENCE) 109 name = "TRAITS AS COVARIATES" 110 CASE (ANALYSE_MULTITRAIT_INCIDENCE_LU) 111 name = "ANALYSE_MULTITRAIT_INCIDENCE_LU" 112 CASE (ANALYSE_MULTITRAIT_INCIDENCE) 113 name = "MULTI-TRAITS ANALYSIS WITH COVARIATE/FIXED EFFECT" 114 CASE (ANALYSE_DEV_1) 115 name = "ANALYSE_DEV_1" 116 CASE (ANALYSE_DEV_2) 117 name = "ANALYSE_DEV_2" 118 CASE DEFAULT 119 call stop_application("Unknown name description for analysis [OPT_CALCUL="//trim(str(opt_calcul))//"]") 120 END SELECT 121 end function get_name_analyse
get_type_data_analyse
[ Top ] [ m_qtlmap_analyse ] [ Subroutines ]
NAME
get_type_data_analyse
DESCRIPTION
Returns the data type which is manageg by the analysis identified by the calcul id
INPUTS
opt_calcul : analysis id
RETURN
Type of data : TYPE_DATA_CONTINUE,TYPE_DATA_COX
SOURCE
134 function get_type_data_analyse(opt_calcul) result (typeData) 135 integer , intent(in) :: opt_calcul 136 integer :: typeData 137 SELECT CASE (opt_calcul) 138 CASE (ANALYSE_UNITRAIT,ANALYSE_UNITRAIT_MODLIN,ANALYSE_MULTITRAIT,ANALYSE_UNITRAIT_LINEAR_LA_HOMO,& 139 ANALYSE_UNITRAIT_LINEAR_LA_HETERO,ANALYSE_UNITRAIT_LDLA,ANALYSE_UNITRAIT_LA,& 140 ANALYSE_UNITRAIT_LD,ANALYSE_UNITRAIT_LDJH, ANALYSE_UNITRAIT_LINEAR_LD_HOMO,ANALYSE_UNITRAIT_LINEAR_LD_HETERO, & 141 ANALYSE_UNITRAIT_LINEAR_LDLA_HETERO,ANALYSE_UNITRAIT_LINEAR_LDLA_HOMO,& 142 ANALYSE_MULTITRAIT_DA,ANALYSE_2QTL_INTERACTION,ANALYSE_UNITRAIT_CONTINGENCE,& 143 ANALYSE_TRAIT_COV_CONTINGENCE,ANALYSE_MULTITRAIT_INCIDENCE,ANALYSE_MULTITRAIT_INCIDENCE_LU,& 144 ANALYSE_DEV_1,ANALYSE_DEV_2) 145 typeData = TYPE_DATA_CONTINUE 146 CASE (ANALYSE_UNITRAIT_MODLIN_COX) 147 typeData = TYPE_DATA_COX 148 CASE DEFAULT 149 call stop_application("Unknown data type description for analysis [OPT_CALCUL="//trim(str(opt_calcul))//"]") 150 END SELECT 151 end function get_type_data_analyse
is_multitrait_analysis
[ Top ] [ m_qtlmap_analyse ] [ Subroutines ]
NAME
is_multitrait_analysis
DESCRIPTION
Returns true if the analysis is a multitrait method, false otherwise
INPUTS
opt_calcul : analysis id
RETURN
boolean : true if the analysis given in parameter is for a set of traits (false otherwise)
SOURCE
164 function is_multitrait_analysis(opt_calcul) result (is_multi) 165 integer , intent(in) :: opt_calcul 166 logical :: is_multi 167 SELECT CASE (opt_calcul) 168 CASE (ANALYSE_UNITRAIT,ANALYSE_UNITRAIT_MODLIN,ANALYSE_UNITRAIT_MODLIN_COX,ANALYSE_UNITRAIT_LINEAR_LA_HOMO,& 169 ANALYSE_UNITRAIT_LINEAR_LA_HETERO,ANALYSE_UNITRAIT_LDLA,ANALYSE_UNITRAIT_LA,ANALYSE_UNITRAIT_LD,& 170 ANALYSE_UNITRAIT_LDJH,ANALYSE_UNITRAIT_LINEAR_LD_HOMO,ANALYSE_UNITRAIT_LINEAR_LD_HETERO, & 171 ANALYSE_UNITRAIT_LINEAR_LDLA_HETERO,ANALYSE_UNITRAIT_LINEAR_LDLA_HOMO,& 172 ANALYSE_2QTL_INTERACTION,ANALYSE_UNITRAIT_CONTINGENCE,& 173 ANALYSE_TRAIT_COV_CONTINGENCE,ANALYSE_DEV_1,ANALYSE_DEV_2) 174 is_multi=.false. 175 return 176 CASE (ANALYSE_MULTITRAIT,ANALYSE_MULTITRAIT_INCIDENCE,ANALYSE_MULTITRAIT_DA,ANALYSE_MULTITRAIT_INCIDENCE_LU) 177 is_multi=.true. 178 return 179 CASE DEFAULT 180 call stop_application("Unknown data type description for analysis [OPT_CALCUL="//trim(str(opt_calcul))//"]") 181 END SELECT 182 183 end function is_multitrait_analysis
need_normalize_data
[ Top ] [ m_qtlmap_analyse ] [ Subroutines ]
NAME
need_normalize_data
DESCRIPTION
Returns true if the analysis need to normalize data
INPUTS
opt_calcul : analysis id
RETURN
boolean : true if the analysis need a data normalization.
SOURCE
197 function need_normalize_data(opt_calcul) result (normalize) 198 integer , intent(in) :: opt_calcul 199 logical :: normalize 200 SELECT CASE (opt_calcul) 201 CASE (ANALYSE_UNITRAIT,ANALYSE_UNITRAIT_MODLIN,ANALYSE_MULTITRAIT,ANALYSE_MULTITRAIT_DA,ANALYSE_UNITRAIT_LINEAR_LA_HOMO,& 202 ANALYSE_UNITRAIT_LINEAR_LA_HETERO,ANALYSE_UNITRAIT_LDLA,ANALYSE_UNITRAIT_LA,ANALYSE_UNITRAIT_LD,& 203 ANALYSE_UNITRAIT_LINEAR_LD_HOMO,ANALYSE_UNITRAIT_LINEAR_LD_HETERO, & 204 ANALYSE_UNITRAIT_LINEAR_LDLA_HETERO,ANALYSE_UNITRAIT_LINEAR_LDLA_HOMO,ANALYSE_UNITRAIT_LDJH,& 205 ANALYSE_2QTL_INTERACTION,ANALYSE_UNITRAIT_CONTINGENCE,& 206 ANALYSE_MULTITRAIT_INCIDENCE,ANALYSE_TRAIT_COV_CONTINGENCE,ANALYSE_MULTITRAIT_INCIDENCE_LU,ANALYSE_DEV_1,& 207 ANALYSE_DEV_2) 208 normalize=.true. 209 return 210 CASE (ANALYSE_UNITRAIT_MODLIN_COX) 211 normalize=.false. 212 return 213 CASE DEFAULT 214 call stop_application("need_normalize_data : Unknown data type description for analysis [OPT_CALCUL="//& 215 trim(str(opt_calcul))//"]") 216 END SELECT 217 218 end function need_normalize_data
opti_multitrait_DA
[ Top ] [ m_qtlmap_analyse ] [ Subroutines ]
NAME
opti_multitrait_DA
DESCRIPTION
LA for a set of traits (without missing data) with a discriminante analysis
INPUTS
opt_sim : contexte of the analysis ( COMMON_ANALYSE | SIMULATION | TRANSCIPTOME_ANALYSE )
OUTPUTS
lrtsol : likelihood ratio test information (see TYPE_LRT_SOLUTION) rhoi : residual correlation matrix finded for each hypothesis
SOURCE
1661 subroutine opti_multitrait_DA(lrtsol,opt_sim) 1662 use m_qtlmap_analyse_gen 1663 use m_qtlmap_analyse_multitrait_DA 1664 type(TYPE_LRT_SOLUTION) , intent(out) :: lrtsol 1665 integer ,intent(in) :: opt_sim 1666 1667 real (kind=dp) , dimension(:,:) ,allocatable :: coeff 1668 real (kind=dp) , dimension(size(carac)) :: coeff_max 1669 real (kind=dp) , dimension(np) :: std_all,ap_all 1670 real (kind=dp) , dimension(nm) :: xmoym1,am1 1671 logical , dimension(nm) :: subphasmm 1672 character(len=LEN_DEF) , dimension(nm) :: submere 1673 logical :: impfem 1674 1675 integer , dimension(np) :: effpp 1676 integer :: i,ic,nmsub,efft,n,ilong,ix,kd,ip,chr 1677 real (kind=dp), dimension (nd) :: yda 1678 integer ,dimension(size(carac)) :: lic 1679 real (kind=dp) :: dxmax 1680 type(TYPE_LRT_SOLUTION) ,dimension(1) :: listLrtSol 1681 type(TYPE_INCIDENCE_SOLUTION) :: incsol 1682 1683 call init_analyse_gen 1684 allocate (coeff(ncar,get_maxnpo())) 1685 1686 if (ncar <= 1) then 1687 call stop_application("can not perform a multitrait analysis on "//trim(str(ncar))//" traits") 1688 end if 1689 1690 do ic=1,ncar 1691 if( natureY(ic) == 'r'.or.natureY(ic) == 'a') cycle 1692 call stop_application("Discriminant multitrait analysis is only available for continuous traits") 1693 end do 1694 1695 call init_analyse_DA_1QTL 1696 1697 coeff=0.d0 1698 1699 call new(1,lrtsol) 1700 ncar=1 1701 lrtsol%lrtmax(0)=-1.d75 1702 dxmax=-1.d6 1703 lrtsol%nxmax=0 1704 lrtsol%chrmax=1 1705 do chr=1,nchr 1706 n=0 1707 ilong=get_ilong(chr) 1708 do ix=0,ilong,pas 1709 n=n+1 1710 call perform_DA_1QTL(chr,n,yda,coeff,get_maxnpo()) 1711 call optinit_da(yda) 1712 call opti_DA_0qtl 1713 call opti_DA_1qtl(chr,ix,n,yda,lrtsol) 1714 if ( lrtsol%chrmax(0) > 0 .and. lrtsol%nxmax(0) > 0 ) then 1715 if(absi(lrtsol%chrmax(0),lrtsol%nxmax(0)) > dxmax) then 1716 dxmax=absi(lrtsol%chrmax(0),lrtsol%nxmax(0)) 1717 coeff_max=coeff(:,n) 1718 end if 1719 call log_mess('at position ix:'//str(ix)//' LRT max:'//str(lrtsol%lrtmax(0))//& 1720 " POS:"//str(absi(lrtsol%chrmax(0),lrtsol%nxmax(0))),DEBUG_DEF) 1721 else 1722 dxmax=0.d0 1723 coeff_max=0.d0 1724 end if 1725 1726 end do 1727 end do 1728 1729 if (opt_sim == COMMON_ANALYSE) then ! Representation graphique 1730 std_all(:)=std 1731 ap_all(:)=ap 1732 call print_start_multitrait_DA 1733 call print_lrt_solution(lrtsol) 1734 call set_solution_hypothesis1(incsol) 1735 call print_incidence_solution(incsol) 1736 call release(incsol) 1737 1738 call print_coeff_linear_combination_max(coeff_max) 1739 call print_paternal_maternal_effect(1,lrtsol) 1740 call print_LRT(lrtsol) 1741 call print_coeff_linear_combination(lrtsol%chrmax(0),get_ilong(lrtsol%chrmax(0)),get_maxnpo(),coeff) 1742 end if 1743 1744 if ( opt_sim /= SIMULATION) then 1745 ! call get_eff_paternal_and_total(effpp,efft) 1746 ! TODO : summary pour multi 1747 end if 1748 1749 ncar=size(carac) 1750 deallocate (coeff) 1751 1752 call end_analyse_DA_1QTL 1753 call end_analyse_gen 1754 1755 end subroutine opti_multitrait_DA
opti_multitrait_modlin
[ Top ] [ m_qtlmap_analyse ] [ Subroutines ]
NAME
opti_multitrait_modlin
DESCRIPTION
LA multi trait with model description (multi qtl version, interaction with fixed effect and qtl are not managed)
INPUTS
nqtl : number of qtl to test. nqtl >= 1 opt_sim : contexte of the analysis ( COMMON_ANALYSE | SIMULATION | TRANSCIPTOME_ANALYSE )
OUTPUTS
lrtsol : likelihood ratio test information (see TYPE_LRT_SOLUTION) listincsol : solution for each hypothesis finded by the analysis (see TYPE_INCIDENCE_SOLUTION)
SOURCE
1027 subroutine opti_multitrait_modlin(lrtsol,listincsol,opt_sim,nqtl,type_resol) 1028 use m_qtlmap_incidence_multi 1029 ! use m_qtlmap_analyse_lin_gen 1030 integer ,intent(in) :: nqtl,type_resol 1031 type(TYPE_LRT_SOLUTION) , intent(inout) ,dimension(nqtl) :: lrtsol 1032 type(TYPE_INCIDENCE_SOLUTION) ,intent(inout) ,dimension(ncar,nqtl+1) :: listincsol 1033 integer ,intent(in) :: opt_sim 1034 1035 integer :: i,j,iq 1036 type(INCIDENCE_GEN_STRUCT) :: workstruct 1037 1038 call init_analyse_gen 1039 call init_workstruct(workstruct,INCIDENCE_TYPE_CLASSIC,OPTIM_HETEROSCEDASTIC,nqtl,& 1040 (opt_sim == COMMON_ANALYSE),(opt_sim == COMMON_ANALYSE),.true.) 1041 if ( type_resol == ANALYSE_MULTITRAIT_INCIDENCE ) then 1042 call opti_0qtl_multi( listincsol(1,1) ,workstruct, nqtl,model_optim_multi_h0,.false.) 1043 ! call opti_0qtl_multi( listincsol(1,1) ,workstruct, nqtl,.false.) 1044 else if ( type_resol == ANALYSE_MULTITRAIT_INCIDENCE_LU ) then 1045 call opti_0qtl_multi( listincsol(1,1) ,workstruct, nqtl,model_optim_multi_h0_LU,.true.) 1046 else 1047 call stop_application("ERROR DEV opti_multitrait_modlin"); 1048 end if 1049 1050 !workstruct%sigsquare(1,:)=workstruct%sigsquareEstime 1051 1052 if (opt_sim == COMMON_ANALYSE) then 1053 call print_start_multitraits() 1054 call print_incidence_solution(listincsol(1,1)) 1055 end if 1056 1057 do i=1,nqtl 1058 workstruct%nqtl=i 1059 if ( type_resol == ANALYSE_MULTITRAIT_INCIDENCE ) then 1060 call opti_nqtl_multi(i,workstruct,listincsol(1,1+i),lrtsol(i),model_optim_multi_hn,nqtl,.false.) 1061 ! call opti_nqtl_multi(i,workstruct,listincsol(1,1+i),lrtsol(i),model_optim_multi_hn,nqtl) 1062 else if ( type_resol == ANALYSE_MULTITRAIT_INCIDENCE_LU ) then 1063 call opti_nqtl_multi(i,workstruct,listincsol(1,1+i),lrtsol(i),model_optim_multi_hn_LU,nqtl,.true.) 1064 else 1065 call stop_application("ERROR DEV opti_multitrait_modlin"); 1066 end if 1067 1068 if (opt_sim == COMMON_ANALYSE) then ! Representation graphique 1069 call print_lrt_solution(lrtsol(i)) 1070 call print_incidence_solution(listincsol(1,1+i)) 1071 1072 !call print_confusion(workstruct%nalert,workstruct%alertQtl,workstruct%corrmax) 1073 ! if ( workstruct(ic)%ntest /= 0 ) then 1074 ! call print_test_nuisances(workstruct(ic)%ntest,workstruct(ic)%listtestnuis) 1075 ! end if 1076 if ( i==1 ) then 1077 !call print_paternal_maternal_effect(ic,lrtsol(ic,1)) 1078 call print_LRT(lrtsol(1)) 1079 end if 1080 1081 if ( i==2 ) then 1082 call create_grid_file_2QTL(lrtsol(2)) 1083 !call print_pat_mat_effect_2QTL(ic,lrtsol(ic,2)) 1084 call print_LRT_2QTL(lrtsol(2)) 1085 end if 1086 1087 end if 1088 end do ! iqtl 1089 1090 call end_analyse_gen 1091 1092 1093 end subroutine opti_multitrait_modlin
opti_multitraits
[ Top ] [ m_qtlmap_analyse ] [ Subroutines ]
NAME
opti_multitraits
DESCRIPTION
LA analysis for a set of traits with a multivariate analysis
INPUTS
opt_sim : contexte of the analysis ( COMMON_ANALYSE | SIMULATION | TRANSCIPTOME_ANALYSE )
OUTPUTS
lrtsol : likelihood ratio test information (see TYPE_LRT_SOLUTION) listincsol : solution for each hypothesis finded by the analysis (see TYPE_INCIDENCE_SOLUTION) rhoi : residual correlation matrix finded for each hypothesis
SOURCE
1592 subroutine opti_multitraits(lrtsol,listincsol,rhoi,opt_sim) 1593 use m_qtlmap_analyse_multitrait 1594 type(TYPE_LRT_SOLUTION) , intent(inout) :: lrtsol 1595 type(TYPE_INCIDENCE_SOLUTION) ,intent(inout) ,dimension(2) :: listincsol 1596 real (kind=dp) ,dimension(size(carac),size(carac)), intent(out) :: rhoi 1597 1598 integer , intent(in) :: opt_sim 1599 1600 logical , dimension(nm) :: subphasm 1601 character(len=LEN_DEF) , dimension(nm) :: submere 1602 logical :: impfem 1603 real (kind=dp) , dimension(nm) :: xmoym1,am1 1604 integer :: nmsub,efft 1605 integer , dimension(np) :: effpp 1606 1607 integer :: ic,i 1608 1609 call init_analyse_gen 1610 1611 if (ncar <= 1) then 1612 call stop_application("can not perform a multitrait analysis on "//trim(str(ncar))//" traits") 1613 end if 1614 1615 do ic=1,ncar 1616 if( natureY(ic) == 'r'.or.natureY(ic) == 'a') cycle 1617 call stop_application("Multitrait analysis is only available for continuous traits") 1618 end do 1619 1620 call init_analyse_multitrait 1621 call optinit_mcar 1622 call opti_mcar_0qtl 1623 call opti_mcar_1qtl(lrtsol) 1624 1625 if (opt_sim == COMMON_ANALYSE .or. opt_sim == TRANSCIPTOME_ANALYSE) then ! Representation graphique 1626 1627 call print_start_multitraits() 1628 call set_solution_hypothesis0(listincsol(1)) 1629 call print_incidence_solution(listincsol(1)) 1630 call print_lrt_solution(lrtsol) 1631 call set_solution_hypothesis1(listincsol(2)) 1632 call print_incidence_solution(listincsol(2)) 1633 1634 call print_paternal_maternal_effect(1,lrtsol) 1635 call print_LRT(lrtsol) 1636 1637 end if 1638 1639 if ( opt_sim /= SIMULATION) then 1640 ! call get_eff_paternal_and_total(effpp,efft) 1641 ! TODO : summary pour multi 1642 end if 1643 1644 call end_analyse_multitrait 1645 call end_analyse_gen 1646 1647 end subroutine opti_multitraits
opti_unitrait
[ Top ] [ m_qtlmap_analyse ] [ Subroutines ]
NAME
opti_unitrait
DESCRIPTION
LA single trait with pre-corrected data (1 or 2 qtl to test)
INPUTS
nbqtl : number of qtl to test. 1 <= nbqtl <= 2 opt_sim : contexte of the analysis ( COMMON_ANALYSE | SIMULATION | TRANSCIPTOME_ANALYSE )
OUTPUTS
lrtsol : likelihood ratio test information (see TYPE_LRT_SOLUTION) listincsol : solution for each hypothesis finded by the analysis (see TYPE_INCIDENCE_SOLUTION)
SOURCE
420 subroutine opti_unitrait(lrtsol,listincsol,opt_sim,nbqtl) 421 use m_qtlmap_analyse_unitrait 422 integer , intent(in) :: nbqtl 423 type(TYPE_LRT_SOLUTION) , intent(inout) ,dimension(ncar,nbqtl) :: lrtsol 424 type(TYPE_INCIDENCE_SOLUTION) , intent(inout) ,dimension(ncar,nbqtl+1) :: listincsol 425 integer , intent(in) :: opt_sim 426 427 integer :: i,ic 428 429 430 431 if ( nbqtl > 2 .or. nbqtl < 1 ) then 432 call stop_application("Devel error : Module UNITRAIT can not perform analysis with nqtl:"//trim(str(nbqtl))) 433 end if 434 435 !$OMP PARALLEL DEFAULT(SHARED) IF (size(carac)>1) 436 call init_analyse_gen 437 call init_analyse_unitrait_1QTL 438 if ( nbqtl == 2) call init_analyse_unitrait_2QTL 439 440 !$OMP DO 441 do ic=1,size(carac) 442 if (.not.(natureY(ic) == 'r'.or.natureY(ic) == 'a')) cycle 443 call log_mess("TRAIT ["//trim(carac(ic))//"]",INFO_DEF) 444 call optinit(ic) 445 call opti_0qtl(ic) 446 call opti_1qtl(ic,lrtsol(ic,1)) 447 448 if ( opt_sim /= SIMULATION) then 449 call set_solution_hypothesis0(ic,listincsol(ic,1)) 450 call set_solution_hypothesis1(ic,listincsol(ic,2)) 451 end if 452 453 if ( nbqtl == 2 ) then 454 call opti_1car_2qtl(ic,lrtsol(ic,2)) 455 if ( opt_sim /= SIMULATION) then 456 call set_solution_hypothesis2(ic,listincsol(ic,3)) 457 end if 458 end if 459 end do 460 !$OMP END DO 461 call end_analyse_unitrait_1QTL 462 if ( nbqtl == 2) call end_analyse_unitrait_2QTL 463 call end_analyse_gen 464 !$OMP END PARALLEL 465 466 if (opt_sim == COMMON_ANALYSE) then ! Representation graphique 467 do ic=1,ncar 468 if (.not.(natureY(ic) == 'r'.or.natureY(ic) == 'a')) cycle 469 call print_start_unitrait(carac(ic)) 470 call print_incidence_solution(listincsol(ic,1)) 471 call print_lrt_solution(lrtsol(ic,1)) 472 call print_incidence_solution(listincsol(ic,2)) 473 call print_paternal_maternal_effect(ic,lrtsol(ic,1)) 474 call print_LRT(lrtsol(ic,1)) 475 if ( nbqtl == 2 ) then 476 call print_lrt_solution(lrtsol(ic,2)) 477 call print_incidence_solution(listincsol(ic,3)) 478 call create_grid_file_2QTL(lrtsol(ic,2)) 479 call print_pat_mat_effect_2QTL(ic,lrtsol(ic,2)) 480 call print_LRT_2QTL(lrtsol(ic,2)) 481 end if 482 end do 483 end if 484 485 486 ! print txt transcriptome format 487 if (opt_sim == TRANSCIPTOME_ANALYSE) then 488 call print_transcriptome_H0(listincsol(:,1)) 489 call print_transcriptome(lrtsol(:,1),listincsol(:,2)) 490 if ( nbqtl == 2) call print_transcriptome(lrtsol(:,2),listincsol(:,3)) 491 end if 492 493 494 if ( opt_sim /= SIMULATION) then 495 call print_summary_analyse(lrtsol(:,1),listincsol(:,2),1) 496 if ( nbqtl == 2) call print_summary_analyse(lrtsol(:,2),listincsol(:,3),2) 497 end if 498 end subroutine opti_unitrait
opti_unitrait_discret
[ Top ] [ m_qtlmap_analyse ] [ Subroutines ]
NAME
opti_unitrait_discret
DESCRIPTION
LA analysis for discrete data with model description (1 qtl)
INPUTS
opt_sim : contexte of the analysis ( COMMON_ANALYSE | SIMULATION | TRANSCIPTOME_ANALYSE )
OUTPUTS
lrtsol : likelihood ratio test information (see TYPE_LRT_SOLUTION) listincsol : solution for each hypothesis finded by the analysis (see TYPE_INCIDENCE_SOLUTION)
SOURCE
1513 subroutine opti_unitrait_discret(lrtsol,listincsol,opt_sim) 1514 use m_qtlmap_analyse_discret_unitrait 1515 use m_qtlmap_analyse_lin_gen 1516 1517 type(TYPE_LRT_SOLUTION) , intent(inout) ,dimension(size(carac),1) :: lrtsol 1518 type(TYPE_INCIDENCE_SOLUTION) ,intent(inout) ,dimension(ncar,2) :: listincsol 1519 1520 integer ,intent(in) :: opt_sim 1521 1522 real (kind=dp) :: fmax 1523 1524 logical , dimension(nm) :: subphasm 1525 1526 integer :: ic,supnbnivest,i 1527 logical :: est_moy = .false. 1528 logical :: est_var = .true. 1529 1530 call init_analyse_gen 1531 do ic=1,size(carac) 1532 if (natureY(ic) /= 'i') cycle 1533 call log_mess("TRAIT ["//trim(carac(ic))//"]",INFO_DEF) 1534 call init_analyse_lin_gen(ic,1) 1535 call init_analyse_discret_unitrait 1536 call optinit(ic) 1537 call opti_0qtl_discret_unitrait(ic) 1538 call opti_1qtl_discret_unitrait(ic,lrtsol(ic,1),fmax,supnbnivest) 1539 1540 if (opt_sim == COMMON_ANALYSE) then ! Representation graphique 1541 call print_start_unitrait(carac(ic)) 1542 ! colinearite donne des indications sur la colinearite entre chacun 1543 ! des effets QTL a estimer et tous les autres effets (avant elimination 1544 ! des effets non estimables) 1545 ! 1546 call confusion('avant',lrtsol(ic,1)%chrmax(0),ic,est_moy,'LA ') 1547 1548 ! Formattage des solutions sous H0 et H1 1549 call set_solution_hypothesis0(ic,listincsol(ic,1)) 1550 call set_solution_hypothesis1(ic,listincsol(ic,2)) 1551 ! Impression du resultat sous H0 1552 call print_incidence_solution(listincsol(ic,1)) 1553 ! Impression de la courbe LRT sous H1 1554 call print_lrt_solution(lrtsol(ic,1)) 1555 ! Impression des resultats sous H1 1556 call print_incidence_solution(listincsol(ic,2)) 1557 ! Effet paternel et maternel sous H1 dans un fichier dedie 1558 call print_paternal_maternal_effect(ic,lrtsol(ic,1)) 1559 ! Affichage de la valeur LRT a chaque position dans un fichier dedie 1560 call print_LRT(lrtsol(ic,1)) 1561 ! 1562 ! test_modlin donne la signification des effets du modele 1563 ! 1564 call test_lin(lrtsol(ic,1)%chrmax(0),ic,est_moy,supnbnivest,fmax,lrtsol(ic,1)%nxmax(0)) 1565 ! 1566 ! confusion indique si'l y a une confusion entre les effets finalement 1567 ! retenus dans le modele 1568 call confusion('apres',lrtsol(ic,1)%chrmax(0),ic,est_moy,'LA ') 1569 end if 1570 1571 call end_analyse_lin_gen 1572 call end_analyse_discret_unitrait 1573 end do 1574 1575 call end_analyse_gen 1576 1577 end subroutine opti_unitrait_discret
opti_unitrait_incidence
[ Top ] [ m_qtlmap_analyse ] [ Subroutines ]
NAME
opti_unitrait_incidence
DESCRIPTION
LA / LD / LDLA single trait with model description with optimised or linearised likelihood.
INPUTS
starticar : start index endicar : end index type_incidence : type of incidence matrix construction (INCIDENCE_TYPE_CLASSIC, INCIDENCE_TYPE_LD, INCIDENCE_TYPE_LDLA, INCIDENCE_TRAIT_COV) type_resolution: type of resolution (LINEAR_HOMOSCEDASTIC, LINEAR_HETEROSCEDASTIC, OPTIM_HETEROSCEDASTIC ) nqtl : number of qtl to test. nqtl >= 1 hdam : LD option, take dam haplotypes in considration opt_sim : contexte of the analysis ( COMMON_ANALYSE | SIMULATION | TRANSCIPTOME_ANALYSE )
OUTPUTS
lrtsol : likelihood ratio test information (see TYPE_LRT_SOLUTION) listincsol : solution for each hypothesis finded by the analysis (see TYPE_INCIDENCE_SOLUTION)
SOURCE
518 ! subroutine opti_unitrait_incidence(lrtsol,listincsol,opt_sim,type_incidence,type_resolution,nqtl,hdam,nopoly) 519 subroutine opti_unitrait_incidence(starticar,endicar,lrtsol,listincsol,opt_sim,type_incidence,& 520 type_resolution,nqtl,hdam,startTraitCov,endTraitCov) 521 use m_qtlmap_incidence_linear 522 use m_qtlmap_incidence_optim 523 use m_qtlmap_analyse_lin_gen 524 integer ,intent(in) :: nqtl,starticar,endicar 525 type(TYPE_LRT_SOLUTION) , intent(inout) ,dimension(ncar,nqtl) :: lrtsol 526 type(TYPE_INCIDENCE_SOLUTION) ,intent(inout) ,dimension(ncar,nqtl+1) :: listincsol 527 ! logical ,intent(in) :: hdam,nopoly 528 logical ,intent(in) :: hdam 529 530 integer ,intent(in) :: opt_sim,type_incidence,type_resolution 531 integer ,optional,intent(in) :: startTraitCov,endTraitCov 532 533 integer :: ic,i,j,ip,jm,stc,etc 534 535 type(INCIDENCE_GEN_STRUCT) ,dimension(:) ,allocatable :: workstruct 536 logical :: islinear 537 538 !start and end trait as cov 539 stc=0 540 etc=0 541 if ( present(startTraitCov) ) stc = startTraitCov 542 if ( present(endTraitCov) ) etc = endTraitCov 543 544 if ( starticar < 1 ) then 545 call stop_application("Dev error [opti_unitrait_incidence] bad def of starticar :"//trim(str(starticar))) 546 end if 547 548 if ( endicar > ncar ) then 549 call stop_application("Dev error [opti_unitrait_incidence] bad def of endicar :"//trim(str(endicar))) 550 end if 551 552 allocate (workstruct(starticar:endicar)) 553 select case (type_resolution) 554 case (LINEAR_HOMOSCEDASTIC) 555 islinear = .true. 556 call init_analyse_linear() 557 case (LINEAR_HETEROSCEDASTIC) 558 islinear = .true. 559 call init_analyse_linear() 560 case (OPTIM_HETEROSCEDASTIC) 561 islinear = .false. 562 case default 563 call stop_application("Devel Error : Unknown incidence analyse type : "//trim(str(type_resolution))) 564 end select 565 566 !$OMP PARALLEL DEFAULT(SHARED) IF (ncar>1) 567 call init_analyse_gen 568 !$OMP DO 569 do ic=starticar,endicar 570 ! if (.not.(natureY(ic) == 'r'.or. natureY(ic) == 'a')) cycle 571 workstruct(ic)%hdam = hdam 572 workstruct(ic)%startTraitAsCov=stc 573 workstruct(ic)%endTraitAsCov=etc 574 ! workstruct(ic)%nopoly = nopoly 575 call log_mess("TRAIT ["//trim(carac(ic))//"]",INFO_DEF) 576 call init_workstruct(workstruct(ic),type_incidence,type_resolution,nqtl,& 577 (opt_sim == COMMON_ANALYSE).and.(type_incidence /= INCIDENCE_TYPE_LD) ,(opt_sim == COMMON_ANALYSE),.false.) 578 579 if (islinear) then 580 call opti_0qtl( ic, listincsol(ic,1),workstruct(ic),model_lin_h0) 581 else 582 call opti_0qtl( ic, listincsol(ic,1),workstruct(ic),model_optim_h0) 583 end if 584 585 do i=1,nqtl 586 workstruct(ic)%nqtl=i 587 if (islinear) then 588 call opti_nqtl(i,ic,workstruct(ic),listincsol(ic,1+i),lrtsol(ic,i),model_lin_hn) 589 else 590 call opti_nqtl(i,ic,workstruct(ic),listincsol(ic,1+i),lrtsol(ic,i),model_optim_hn) 591 end if 592 end do ! iqtl 593 end do ! ic 594 !$OMP END DO 595 call end_analyse_gen 596 !$OMP END PARALLEL 597 598 599 if (opt_sim == COMMON_ANALYSE) then ! Representation graphique 600 do ic=starticar,endicar 601 call print_start_unitrait(carac(ic)) 602 call print_incidence_solution(listincsol(ic,1)) 603 do i=1,nqtl 604 call print_lrt_solution(lrtsol(ic,i)) 605 call print_incidence_solution(listincsol(ic,1+i)) 606 call print_confusion(workstruct(ic)%nalert,workstruct(ic)%alertQtl,workstruct(ic)%corrmax) 607 if ( workstruct(ic)%ntest /= 0 ) then 608 call print_test_nuisances(workstruct(ic)%ntest,workstruct(ic)%listtestnuis) 609 end if 610 if ( i==1 ) then 611 call print_paternal_maternal_effect(ic,lrtsol(ic,1)) 612 call print_LRT(lrtsol(ic,1)) 613 end if 614 if ( i==2 ) then 615 call create_grid_file_2QTL(lrtsol(ic,2)) 616 call print_pat_mat_effect_2QTL(ic,lrtsol(ic,2)) 617 call print_LRT_2QTL(lrtsol(ic,2)) 618 end if 619 end do ! iqtl 620 end do 621 end if 622 623 if ( opt_sim /= SIMULATION .and. (type_incidence /= INCIDENCE_TYPE_LD) ) then 624 625 do i=1,nqtl 626 call print_summary_analyse(lrtsol(:,i),listincsol(:,i+1),i,starticar,endicar) 627 end do 628 end if 629 630 ! print txt transcriptome format 631 if (opt_sim == TRANSCIPTOME_ANALYSE) then 632 call print_transcriptome_H0(listincsol(:,1)) 633 do i=1,nqtl 634 call print_transcriptome(lrtsol(:,i),listincsol(:,i+1)) 635 end do 636 end if 637 638 select case (type_resolution) 639 case (LINEAR_HOMOSCEDASTIC) 640 call end_analyse_linear 641 case (LINEAR_HETEROSCEDASTIC) 642 call end_analyse_linear 643 case (OPTIM_HETEROSCEDASTIC) 644 islinear = .false. 645 case default 646 call stop_application("Devel Error : Unknown incidence analyse type : "//trim(str(type_resolution))) 647 end select 648 649 do ic=starticar,endicar 650 if (.not.(natureY(ic) == 'r'.or.natureY(ic) == 'a')) cycle 651 call release_ws(workstruct(ic)) 652 end do 653 654 deallocate (workstruct) 655 656 end subroutine opti_unitrait_incidence
opti_unitrait_incidence_cuda
[ Top ] [ m_qtlmap_analyse ] [ Subroutines ]
NAME
opti_unitrait_incidence_cuda
DESCRIPTION
LA / LD / LDLA single trait with model description with optimised or linearised likelihood.
INPUTS
type_incidence : type of incidence matrix construction (INCIDENCE_TYPE_CLASSIC, INCIDENCE_TYPE_LD, INCIDENCE_TYPE_LDLA, INCIDENCE_TRAIT_COV) type_resolution: type of resolution (LINEAR_HOMOSCEDASTIC, LINEAR_HETEROSCEDASTIC, OPTIM_HETEROSCEDASTIC ) nqtl : number of qtl to test. nqtl >= 1 hdam : LD option, take dam haplotypes in considration opt_sim : contexte of the analysis ( COMMON_ANALYSE | SIMULATION | TRANSCIPTOME_ANALYSE )
OUTPUTS
lrtsol : likelihood ratio test information (see TYPE_LRT_SOLUTION) listincsol : solution for each hypothesis finded by the analysis (see TYPE_INCIDENCE_SOLUTION)
SOURCE
895 #ifdef MANAGE_CUDA 896 subroutine opti_unitrait_incidence_cuda(nsim,lrtsol,listincsol,listLrtSolSim,opt_sim,type_incidence,& 897 type_resolution,nqtl,ySIMUL,hdam) 898 use m_qtlmap_incidence_linear 899 use m_qtlmap_incidence_optim 900 use m_qtlmap_analyse_lin_gen 901 integer ,intent(in) :: nqtl,opt_sim,nsim 902 type(TYPE_LRT_SOLUTION) , intent(inout) ,dimension(ncar,nqtl) :: lrtsol 903 type(TYPE_INCIDENCE_SOLUTION) ,intent(inout) ,dimension(ncar,nqtl+1) :: listincsol 904 type(TYPE_LRT_SOLUTION) , intent(inout) ,dimension(nsim,ncar,nqtl) :: listLrtSolSim 905 real (kind=dp), dimension (nsim,ncar,nd) , intent(in) :: ySIMUL 906 logical ,intent(in) :: hdam 907 908 integer ,intent(in) :: type_incidence,type_resolution 909 910 integer :: ic,i,j,ip,jm,kd 911 real (kind=dp), dimension (:,:) ,allocatable :: AllsigQtlMoinUn 912 913 type(INCIDENCE_GEN_STRUCT) ,dimension(:) ,allocatable :: workstruct 914 logical :: islinear 915 integer :: nkd,sizeF(np) 916 real (kind=dp), dimension (:,:),allocatable :: YCOMPACT 917 918 allocate (YCOMPACT(nsim+1,nd)) 919 allocate (workstruct(ncar)) 920 allocate (AllsigQtlMoinUn(nsim+1,np)) 921 922 call init_analyse_linear() 923 call init_analyse_gen 924 do ic=1,size(carac) 925 YCOMPACT=0.d0 926 AllsigQtlMoinUn=0.d0 927 nkd=0 928 do ip=1,np 929 do jm=nmp(ip)+1,nmp(ip+1) 930 do kd=ndm(jm)+1,ndm(jm+1) 931 if (presentc(ic,kd)) then 932 nkd=nkd+1 933 YCOMPACT(1,nkd) = y(ic,kd) 934 YCOMPACT(2:nsim+1,nkd)= ySIMUL(:,ic,kd) 935 !YCOMPACT(2:nsim,nkd)= y(ic,kd) ! pour les test 936 end if 937 end do !kd 938 end do ! jm 939 end do ! ip 940 941 workstruct(ic)%hdam = hdam 942 call log_mess("TRAIT ["//trim(carac(ic))//"]",INFO_DEF) 943 call init_workstruct(workstruct(ic),type_incidence,type_resolution,nqtl,& 944 (type_incidence /= INCIDENCE_TYPE_LD) ,.true.,.false.) 945 946 call opti_0qtl_cuda( nsim,ic, listincsol(ic,1),workstruct(ic),YCOMPACT,sizeF,AllsigQtlMoinUn,model_lin_h0) 947 !call opti_0qtl( ic, listincsol(ic,1),workstruct(ic),model_lin_h0) 948 949 do i=1,nqtl 950 workstruct(ic)%nqtl=i 951 call opti_nqtl_cuda(nsim,i,ic,workstruct(ic),listincsol(ic,1+i),lrtsol(ic,i),& 952 YCOMPACT,listLrtSolSim,sizeF,AllsigQtlMoinUn,model_lin_hn) 953 954 ! call opti_nqtl(i,ic,workstruct(ic),listincsol(ic,1+i),lrtsol(ic,i),model_lin_hn) 955 end do ! iqtl 956 end do ! ic 957 call end_analyse_gen 958 959 call release_allocated_internal_sig 960 961 if (opt_sim /= TRANSCIPTOME_ANALYSE) then ! Representation graphique 962 do ic=1,ncar 963 call print_start_unitrait(carac(ic)) 964 call print_incidence_solution(listincsol(ic,1)) 965 do i=1,nqtl 966 call print_lrt_solution(lrtsol(ic,i)) 967 call print_incidence_solution(listincsol(ic,1+i)) 968 call print_confusion(workstruct(ic)%nalert,workstruct(ic)%alertQtl,workstruct(ic)%corrmax) 969 if ( workstruct(ic)%ntest /= 0 ) then 970 call print_test_nuisances(workstruct(ic)%ntest,workstruct(ic)%listtestnuis) 971 end if 972 if ( i==1 ) then 973 call print_paternal_maternal_effect(ic,lrtsol(ic,1)) 974 call print_LRT(lrtsol(ic,1)) 975 end if 976 if ( i==2 ) then 977 call create_grid_file_2QTL(lrtsol(ic,2)) 978 call print_pat_mat_effect_2QTL(ic,lrtsol(ic,2)) 979 call print_LRT_2QTL(lrtsol(ic,2)) 980 end if 981 end do ! iqtl 982 end do 983 end if 984 985 if ( (type_incidence /= INCIDENCE_TYPE_LD) ) then 986 do i=1,nqtl 987 call print_summary_analyse(lrtsol(:,i),listincsol(:,i+1),i) 988 end do 989 end if 990 991 ! print txt transcriptome format 992 if (opt_sim == TRANSCIPTOME_ANALYSE) then 993 call print_transcriptome_H0(listincsol(:,1)) 994 do i=1,nqtl 995 call print_transcriptome(lrtsol(:,i),listincsol(:,i+1)) 996 end do 997 end if 998 999 call end_analyse_linear 1000 1001 do ic=1,size(carac) 1002 if (.not.(natureY(ic) == 'r'.or.natureY(ic) == 'a')) cycle 1003 call release_ws(workstruct(ic)) 1004 end do 1005 1006 deallocate (workstruct) 1007 deallocate (YCOMPACT) 1008 deallocate (AllsigQtlMoinUn) 1009 1010 end subroutine opti_unitrait_incidence_cuda 1011 #endif
opti_unitrait_interaction
[ Top ] [ m_qtlmap_analyse ] [ Subroutines ]
NAME
opti_unitrait_interaction
DESCRIPTION
*** DEV interaction qtl analysis ****
SOURCE
1102 subroutine opti_unitrait_interaction(lrtsol,opt_sim) 1103 use m_qtlmap_incidence_linear 1104 use m_qtlmap_incidence_optim 1105 use m_qtlmap_analyse_lin_gen 1106 type(TYPE_LRT_SOLUTION) , intent(inout) ,dimension(size(carac),2) :: lrtsol 1107 integer ,intent(in) :: opt_sim 1108 1109 integer :: type_resolution = LINEAR_HOMOSCEDASTIC 1110 integer :: type_incidence = INCIDENCE_INTERACTION 1111 !type(TYPE_INCIDENCE_SOLUTION) ,intent(inout) ,dimension(ncar,3) :: listincsol 1112 1113 integer :: ic,supnbnivest,typeModel,i,j,ip 1114 1115 real (kind=dp) ,dimension(0:2,np) :: sigsquare 1116 real (kind=dp) ,dimension(np) :: sigsquareEstime 1117 1118 !1 ->H0, 2->H1, 3->H2, 4->H3 (avec interaction) 1119 type(TYPE_INCIDENCE_SOLUTION) ,dimension(size(carac),4) :: listincsol 1120 type(TYPE_LRT_SOLUTION) ,dimension(size(carac)) :: lrtsolInteraction 1121 1122 type(INCIDENCE_GEN_STRUCT) ,dimension(:) ,allocatable :: workstruct 1123 logical :: islinear 1124 1125 1126 call init_analyse_gen 1127 allocate (workstruct(ncar)) 1128 1129 select case (type_resolution) 1130 case (LINEAR_HOMOSCEDASTIC) 1131 islinear = .true. 1132 call init_analyse_linear() 1133 case (LINEAR_HETEROSCEDASTIC) 1134 islinear = .true. 1135 call init_analyse_linear() 1136 case (OPTIM_HETEROSCEDASTIC) 1137 islinear = .false. 1138 case default 1139 call stop_application("Devel Error : Unknown incidence analyse type : "//trim(str(type_resolution))) 1140 end select 1141 1142 do ic=1,size(carac) 1143 if (.not.(natureY(ic) == 'r'.or. natureY(ic) == 'a')) cycle 1144 sigsquare=0.d0 1145 1146 call init_workstruct(workstruct(ic),type_incidence,type_resolution,2,& 1147 (opt_sim == COMMON_ANALYSE),(opt_sim == COMMON_ANALYSE),.false.) 1148 if (islinear) then 1149 call opti_0qtl( ic, listincsol(ic,1),workstruct(ic),model_lin_h0) 1150 else 1151 call opti_0qtl( ic, listincsol(ic,1),workstruct(ic),model_optim_h0) 1152 end if 1153 1154 if (opt_sim == COMMON_ANALYSE) then 1155 call print_start_unitrait(carac(ic)) 1156 call print_incidence_solution(listincsol(ic,1)) 1157 end if 1158 1159 do i=1,2 1160 workstruct(ic)%nqtl=i 1161 if (islinear) then 1162 call opti_nqtl(i,ic,workstruct(ic),listincsol(ic,1+i),lrtsol(ic,i),model_lin_hn) 1163 else 1164 call opti_nqtl(i,ic,workstruct(ic),listincsol(ic,1+i),lrtsol(ic,i),model_optim_hn) 1165 end if 1166 1167 if (opt_sim == COMMON_ANALYSE) then ! Representation graphique 1168 call print_lrt_solution(lrtsol(ic,i)) 1169 call print_incidence_solution(listincsol(ic,1+i)) 1170 call print_confusion(workstruct(ic)%nalert,workstruct(ic)%alertQtl,workstruct(ic)%corrmax) 1171 if ( workstruct(ic)%ntest /= 0 ) then 1172 call print_test_nuisances(workstruct(ic)%ntest,workstruct(ic)%listtestnuis) 1173 end if 1174 if ( i==1 ) then 1175 call print_paternal_maternal_effect(ic,lrtsol(ic,1)) 1176 call print_LRT(lrtsol(ic,1)) 1177 end if 1178 1179 if ( i==2 ) then 1180 call create_grid_file_2QTL(lrtsol(ic,2)) 1181 call print_pat_mat_effect_2QTL(ic,lrtsol(ic,2)) 1182 call print_LRT_2QTL(lrtsol(ic,2)) 1183 end if 1184 1185 end if 1186 end do ! i 1187 1188 !Test avec interaction seulement lineaire la resolution.....a voir par la suite... 1189 1190 call opti_2qtl_linear_interaction(ic,workstruct(ic),listincsol(ic,4),lrtsolInteraction(ic),typeModel) 1191 1192 if (opt_sim == COMMON_ANALYSE) then ! Representation graphique 1193 call print_lrt_solution(lrtsol(ic,2)) 1194 call print_incidence_solution(listincsol(ic,4)) 1195 call create_grid_file_2QTL(lrtsol(ic,2)) 1196 call print_pat_mat_effect_2QTL(ic,lrtsol(ic,2)) 1197 call print_LRT_2QTL(lrtsol(ic,2)) 1198 end if 1199 1200 ! end if 1201 end do ! ic 1202 call end_analyse_linear 1203 1204 do ic=1,size(carac) 1205 if (.not.(natureY(ic) == 'r'.or.natureY(ic) == 'a')) cycle 1206 do i=1,3 1207 call release(listincsol(ic,i)) 1208 end do 1209 end do 1210 1211 do ic=1,size(carac) 1212 if (.not.(natureY(ic) == 'r'.or.natureY(ic) == 'a')) cycle 1213 call release_ws(workstruct(ic)) 1214 end do 1215 1216 deallocate (workstruct) 1217 call end_analyse_gen 1218 1219 1220 end subroutine opti_unitrait_interaction
opti_unitrait_modlin
[ Top ] [ m_qtlmap_analyse ] [ Subroutines ]
NAME
opti_unitrait_modlin
DESCRIPTION
LA single trait with model description (1 qtl)
INPUTS
opt_sim : contexte of the analysis ( COMMON_ANALYSE | SIMULATION | TRANSCIPTOME_ANALYSE )
OUTPUTS
lrtsol : likelihood ratio test information (see TYPE_LRT_SOLUTION) listincsol : solution for each hypothesis finded by the analysis (see TYPE_INCIDENCE_SOLUTION)
SOURCE
1234 subroutine opti_unitrait_modlin(lrtsol,listincsol,opt_sim) 1235 use m_qtlmap_analyse_modlin 1236 !#ifdef MANAGE_OMP 1237 ! use m_qtlmap_analyse_modlin_omp 1238 !#endif 1239 use m_qtlmap_analyse_lin_gen 1240 1241 type(TYPE_LRT_SOLUTION) , intent(inout) ,dimension(size(carac),1) :: lrtsol 1242 integer ,intent(in) :: opt_sim 1243 type(TYPE_INCIDENCE_SOLUTION) ,intent(inout) ,dimension(ncar,2) :: listincsol 1244 real (kind=dp) :: fmax 1245 1246 integer :: ic,supnbnivest 1247 !modif carole : l'estimationd e la moyenne est parametrable 1248 logical :: est_moy = .true. 1249 logical :: est_var = .true. 1250 integer :: i 1251 1252 !$OMP PARALLEL DEFAULT(SHARED) PRIVATE(fmax,supnbnivest) IF ( opt_sim /= COMMON_ANALYSE .and. (size(carac)>1) ) 1253 call init_analyse_gen 1254 call init_analyse_modlin 1255 1256 !$OMP DO 1257 do ic=1,size(carac) 1258 if (.not.(natureY(ic) == 'r'.or.natureY(ic) == 'a')) cycle 1259 call log_mess("TRAIT ["//trim(carac(ic))//"]",INFO_DEF) 1260 call init_analyse_lin_gen(ic,1) 1261 call optinit(ic) 1262 call opti_0qtl_modlin(ic) 1263 call opti_1qtl_modlin(ic,lrtsol(ic,1),fmax,supnbnivest) 1264 1265 1266 if ( opt_sim /= SIMULATION) then 1267 ! Formattage des solutions sous H0 et H1 1268 call set_solution_hypothesis0(ic,listincsol(ic,1)) 1269 1270 call set_solution_hypothesis1(ic,listincsol(ic,2)) 1271 end if 1272 1273 !a cause de test_lin ,on ne peut pas sortir l affichage de la boucle de traitement 1274 !cette region devient donc critique pour afficher correctement et succesivement par chaque thread les resultats.... 1275 if (opt_sim == COMMON_ANALYSE) then ! Representation graphique 1276 !$OMP CRITICAL 1277 call print_start_unitrait(carac(ic)) 1278 ! colinearite donne des indications sur la colinearite entre chacun 1279 ! des effets QTL a estimer et tous les autres effets (avant elimination 1280 ! des effets non estimables) 1281 ! 1282 call confusion('avant',lrtsol(ic,1)%chrmax(0),ic,est_moy,'LA ') 1283 ! Impression du resultat sous H0 1284 call print_incidence_solution(listincsol(ic,1)) 1285 ! Impression de la courbe LRT sous H1 1286 call print_lrt_solution(lrtsol(ic,1)) 1287 ! Impression des resultats sous H1 1288 call print_incidence_solution(listincsol(ic,2)) 1289 ! Effet paternel et maternel sous H1 dans un fichier dedie 1290 call print_paternal_maternal_effect(ic,lrtsol(ic,1)) 1291 ! Affichage de la valeur LRT a chaque position dans un fichier dedie 1292 call print_LRT(lrtsol(ic,1)) 1293 ! 1294 ! test_modlin donne la signification des effets du modele 1295 ! 1296 call test_lin(lrtsol(ic,1)%chrmax(0),ic,est_moy,supnbnivest,fmax,lrtsol(ic,1)%nxmax(0)) 1297 ! 1298 ! confusion indique si'l y a une confusion entre les effets finalement 1299 ! retenus dans le modele 1300 call confusion('apres',lrtsol(ic,1)%chrmax(0),ic,est_moy,'LA ') 1301 !$OMP END CRITICAL 1302 end if 1303 call end_analyse_lin_gen 1304 end do 1305 !$OMP END DO 1306 call end_analyse_modlin 1307 call end_analyse_gen 1308 !$OMP END PARALLEL 1309 1310 if ( opt_sim /= SIMULATION) then 1311 call print_summary_analyse(lrtsol(:,1),listincsol(:,2),1) 1312 end if 1313 1314 ! print txt transcriptome format 1315 if (opt_sim == TRANSCIPTOME_ANALYSE) then 1316 call print_transcriptome_H0(listincsol(:,1)) 1317 call print_transcriptome(lrtsol(:,1),listincsol(:,2)) 1318 end if 1319 1320 end subroutine opti_unitrait_modlin
opti_unitrait_modlin_cox
[ Top ] [ m_qtlmap_analyse ] [ Subroutines ]
NAME
opti_unitrait_modlin_cox
DESCRIPTION
LA analysis with cox model
INPUTS
opt_sim : contexte of the analysis ( COMMON_ANALYSE | SIMULATION | TRANSCIPTOME_ANALYSE )
OUTPUTS
lrtsol : likelihood ratio test information (see TYPE_LRT_SOLUTION) listincsol : solution for each hypothesis finded by the analysis (see TYPE_INCIDENCE_SOLUTION)
SOURCE
1428 subroutine opti_unitrait_modlin_cox(lrtsol,listincsol,opt_sim) 1429 use m_qtlmap_analyse_modlin_cox 1430 use m_qtlmap_analyse_lin_gen 1431 1432 type(TYPE_LRT_SOLUTION) , intent(inout) ,dimension(size(carac),1) :: lrtsol 1433 type(TYPE_INCIDENCE_SOLUTION) ,intent(inout) ,dimension(ncar,2) :: listincsol 1434 integer ,intent(in) :: opt_sim 1435 1436 real (kind=dp) :: fmax 1437 integer :: ic,supnbnivest,npar,i 1438 logical :: est_moy = .false. 1439 logical :: est_var = .false. 1440 1441 call init_analyse_gen 1442 do ic=1,size(carac) 1443 if (natureY(ic) /= 'r') cycle 1444 call log_mess("TRAIT ["//trim(carac(ic))//"]",INFO_DEF) 1445 call init_analyse_lin_gen(ic,1) 1446 !npar is modified because it is equal to 3*np + 2nm in m_qtlmap_analyse_lin_gen.f95 1447 npar=np+nm 1448 call init_analyse_modlin_cox 1449 call optinit(ic) 1450 call calcul_rang(ic) 1451 call opti_0qtl_modlin_cox(ic,est_moy) 1452 call opti_1qtl_modlin_cox(ic,lrtsol(ic,1),fmax,supnbnivest) 1453 !call stop_application("Survival Analysis with a Cox model (option=5) is in development ") 1454 1455 if (opt_sim == COMMON_ANALYSE) then ! Representation graphique 1456 call print_start_unitrait(carac(ic)) 1457 ! colinearite donne des indications sur la colinearite entre chacun 1458 ! des effets QTL a estimer et tous les autres effets (avant elimination 1459 ! des effets non estimables) 1460 call confusion('avant',lrtsol(ic,1)%chrmax(0),ic,.false.,'LA ') 1461 ! Impression de la courbe LRT sous H1 1462 call print_lrt_solution(lrtsol(ic,1)) 1463 1464 ! Formattage des solutions sous H0 et H1 1465 call set_solution_hypothesis0(ic,listincsol(ic,1)) 1466 call set_solution_hypothesis1(ic,listincsol(ic,2)) 1467 ! Impression du resultat sous H0 1468 call print_incidence_solution_risk_factor(listincsol(ic,1)) 1469 ! Impression des resultats sous H1 1470 call print_incidence_solution_risk_factor(listincsol(ic,2)) 1471 ! Effet paternel et maternel sous H1 dans un fichier dedie 1472 call print_paternal_maternal_effect(ic,lrtsol(ic,1)) 1473 ! Affichage de la valeur LRT a chaque position dans un fichier dedie 1474 call print_LRT(lrtsol(ic,1)) 1475 ! 1476 ! test_modlin donne la signification des effets du modele 1477 ! 1478 call test_lin_cox(lrtsol(ic,1)%chrmax(0),ic,est_moy,supnbnivest,fmax,lrtsol(ic,1)%nxmax(0)) 1479 ! 1480 ! confusion indique si'l y a une confusion entre les effets finalement 1481 ! retenus dans le modele 1482 call confusion('apres',lrtsol(ic,1)%chrmax(0),ic,est_moy,'LA ') 1483 end if 1484 call end_analyse_lin_gen 1485 call end_analyse_modlin_cox 1486 enddo 1487 1488 do ic=1,size(carac) 1489 if (natureY(ic) /= 'r') cycle 1490 do i=1,2 1491 call release(listincsol(ic,i)) 1492 end do 1493 end do 1494 1495 !TODO: Olivier il faut mettre le call print_summary_analyse(efft,effpp,lrtsol(:,1),listincsol(:,2),1) je n'y arrive pas!!!!! (carole) 1496 1497 call end_analyse_gen 1498 end subroutine opti_unitrait_modlin_cox
opti_unitrait_modlin_new
[ Top ] [ m_qtlmap_analyse ] [ Subroutines ]
NAME
opti_unitrait_modlin_new
DESCRIPTION
LA / LD / LDLA /LDJH single trait with model description (1 qtl) optimisation method
INPUTS
option_var : "hete" : heteroscedastic, "homo" : homoscedastic option_anal : "LDLA", "LD ", "LA ","LDJH" hdam : LD option, take dam haplotypes in considration opt_sim : contexte of the analysis ( COMMON_ANALYSE | SIMULATION | TRANSCIPTOME_ANALYSE )
OUTPUTS
lrtsol : likelihood ratio test information (see TYPE_LRT_SOLUTION) listincsol : solution for each hypothesis finded by the analysis (see TYPE_INCIDENCE_SOLUTION)
SOURCE
1337 subroutine opti_unitrait_modlin_new(lrtsol,listincsol,option_var,option_anal,hdam,opt_sim) 1338 use m_qtlmap_analyse_modlin_ldla 1339 use m_qtlmap_analyse_lin_gen 1340 1341 type(TYPE_LRT_SOLUTION) , intent(inout) ,dimension(size(carac),1) :: lrtsol 1342 type(TYPE_INCIDENCE_SOLUTION) ,intent(inout) ,dimension(ncar,2) :: listincsol 1343 character(len=4) , intent(in) :: option_var,option_anal 1344 integer ,intent(in) :: opt_sim 1345 logical ,intent(in) :: hdam 1346 1347 real (kind=dp) , dimension(size(carac),np) :: std_all,ap_all 1348 logical , dimension(nm) :: subphasm 1349 real(kind=dp) :: fmax 1350 integer :: ic,supnbnivest,nposx,CHR 1351 !modif carole : l'estimationd e la moyenne est parametrable 1352 logical :: est_moy = .true. 1353 logical :: est_var = .true. 1354 logical :: qtl,hsire 1355 1356 hsire=( (option_anal == 'LD ') .or. (option_anal == 'LDLA') .or. (option_anal == 'LDJH') ) 1357 qtl = ( (option_anal == 'LA ') .or. (option_anal == 'LDLA') .or. (option_anal == 'LDJH') ) 1358 1359 call init_analyse_gen 1360 call init_analyse_modlin_ldla 1361 1362 do ic=1,size(carac) 1363 if (.not.(natureY(ic) == 'r'.or.natureY(ic) == 'a')) cycle 1364 call log_mess("TRAIT ["//trim(carac(ic))//"]",INFO_DEF) 1365 call init_analyse_lin_gen(ic,1) 1366 call optinit(ic) 1367 call opti_0qtl_modlin_ldla(ic,.true.,option_var) 1368 call opti_1qtl_modlin_ldla(ic,lrtsol(ic,1),fmax,supnbnivest,.true.,option_var,option_anal,hdam) 1369 1370 if ( opt_sim /= SIMULATION) then 1371 ! Formattage des solutions sous H0 et H1 1372 call set_solution_hypothesis0(ic,(option_var=="hete"),listincsol(ic,1)) 1373 call set_solution_hypothesis1(ic,(option_var=="hete"),qtl,hsire,hdam,listincsol(ic,2)) 1374 end if 1375 1376 if (opt_sim == COMMON_ANALYSE) then ! Representation graphique 1377 call print_start_unitrait(carac(ic)) 1378 ! colinearite donne des indications sur la colinearite entre chacun 1379 ! des effets QTL a estimer et tous les autres effets (avant elimination 1380 ! des effets non estimables) 1381 ! 1382 call confusion('avant',lrtsol(ic,1)%chrmax(0),ic,est_moy,option_anal) 1383 ! Impression du resultat sous H0 1384 call print_incidence_solution(listincsol(ic,1)) 1385 ! Impression de la courbe LRT sous H1 1386 call print_lrt_solution(lrtsol(ic,1)) 1387 ! Impression des resultats sous H1 1388 call print_incidence_solution(listincsol(ic,2)) 1389 ! Effet paternel et maternel sous H1 dans un fichier dedie 1390 call print_paternal_maternal_effect(ic,lrtsol(ic,1)) 1391 ! Affichage de la valeur LRT a chaque position dans un fichier dedie 1392 call print_LRT(lrtsol(ic,1)) 1393 ! 1394 ! test_modlin donne la signification des effets du modele 1395 ! 1396 ! call test_lin(lrtsol(ic,1)%chrmax(0),ic,est_moy,supnbnivest,fmax,lrtsol(ic,1)%nxmax(0)) 1397 ! 1398 ! confusion indique si'l y a une confusion entre les effets finalement 1399 ! retenus dans le modele 1400 call confusion('apres',lrtsol(ic,1)%chrmax(0),ic,est_moy,'LA ') 1401 end if 1402 1403 call end_analyse_lin_gen 1404 end do 1405 1406 ! if ( opt_sim /= OPT_SIMULATION) then 1407 ! call get_eff_paternal_and_total(effp,efft) 1408 ! call print_summary_0versus1(size(carac),carac,sigt,size(pere),pere,efft,effp,std_all,ap_all,xlrmax0_1,dxma0_1) 1409 ! end if 1410 1411 call end_analyse_modlin_ldla 1412 call end_analyse_gen 1413 1414 end subroutine opti_unitrait_modlin_new
quantile
[ Top ] [ m_qtlmap_analyse ] [ Subroutines ]
NAME
quantile
DESCRIPTION
computes quantiles (Harrell & Davis 1982) and write a result file.
INPUTS
simula_file : name of file to print the output nsim : number of simulation opt_calcul : analysis id opt_qtl : number of qtl hypothesis (analysis can be not compatible with this option) nqtl : number of qtl simulated lrtsol : likelihood ratio test information (see TYPE_LRT_SOLUTION)
NOTES
use : MATH_QTLMAP_M01CAF, MATH_QTLMAP_G01AAF, QUANTILE1, print_resume_simulation, print_resume_simulation_2
SOURCE
1773 subroutine quantile(simula_file,nsim,opt_calcul,opt_qtl,nqtl,lrtsol) 1774 character(len=LENGTH_MAX_FILE) , intent(in) :: simula_file 1775 integer , intent(in) :: nsim,opt_calcul,opt_qtl,nqtl ! nqtl => le nombre de qtl definit dans la simulation 1776 ! opt_qtl => nombre de qtl a detecter pour le calcul 1777 type(TYPE_LRT_SOLUTION) , dimension(nsim,ncar,opt_qtl),intent(in) :: lrtsol 1778 1779 real (kind=dp) :: prob(16) 1780 real (kind=dp) :: x(nsim),wt(nsim) 1781 real (kind=dp) ,dimension (ncar,16) :: z 1782 real (kind=dp) ,dimension (1,16) :: z2 1783 1784 CHARACTER(len=LEN_BUFFER_WORD) :: temp,nameT,carac_t(1) 1785 integer :: nc,ifail,i,im,j,ii,iqtl,iwt 1786 real (kind=dp) :: xmin1,xmax1,ymu1,sig1,s21,s31,pg,wtsum,p_t,do 1787 1788 data prob /.1d0,.2d0,.3d0,.4d0,.5d0,.6d0,.7d0,.8d0,.9d0,.95d0,& 1789 .99d0,.995d0,.9973d0,.999d0,.99947d0,.999947d0/ 1790 1791 if ( nsim < 2 ) return 1792 1793 nc=size(carac) 1794 nameT="" 1795 if ( is_multitrait_analysis(opt_calcul) ) then 1796 nc = 1 1797 nameT="* All traits *" 1798 end if 1799 1800 ifail=0 1801 iwt = 0 1802 ! pour tous les lrts trouve sous avec qtl==opt_qtl 1803 do iqtl=1,opt_qtl ! ? / opt_qtl 1804 do im=0,iqtl-1 ! on ecrit tous les test H: im / opt_qtl 1805 temp=trim(str(im))//'vs'//trim(str(iqtl))//'Q' 1806 do i=1,nc 1807 do j=1,nsim 1808 x(j)=lrtsol(j,i,iqtl)%lrtmax(im) 1809 end do 1810 if (nameT=="") nameT=trim(carac(i)) 1811 1812 !Sort X from index 1 to ns ('a' ascendant).... 1813 call MATH_QTLMAP_M01CAF(X,1,nsim,'A',IFAIL) 1814 !Computing Mean, variance, skewness, kurtosis Minimum Maximum 1815 call MATH_QTLMAP_G01AAF(nsim,x,IWT,WT,YMU1,SIG1,S21,S31,XMIN1,XMAX1,WTSUM,IFAIL) 1816 do ii=9,16 1817 z(i,ii)=QUANTILE1(nsim,prob(ii),X) 1818 end do 1819 1820 call print_resume_simulation(nameT,temp,nsim,& 1821 YMU1,SIG1,S21,S31,XMIN1,XMAX1,prob,z(i,:),9,16) 1822 nameT="" 1823 end do 1824 if ( .not. is_multitrait_analysis(opt_calcul) ) then 1825 call print_resume_simulation_2(temp,ncar,carac,z) 1826 else 1827 carac_t(1)=trim(nameT) 1828 z2(1,:)=z(1,:) 1829 call print_resume_simulation_2(temp,1,carac_t,z2) 1830 end if 1831 end do 1832 end do 1833 1834 end subroutine quantile
QUANTILE1
[ Top ] [ m_qtlmap_analyse ] [ Subroutines ]
NAME
QUANTILE1
DESCRIPTION
CALCUL DES QUANTILES D'UNE DISTRIBUTION ORDONNEE
INPUTS
n : p : x :
RETURN
NOTES
use : W
SOURCE
1851 FUNCTION QUANTILE1(n,p,X) result (Q) 1852 integer , intent(in) :: n 1853 real (kind=dp) ,intent(in) :: p 1854 real (kind=dp) , dimension(:),intent(in) :: X 1855 1856 real (kind=dp) :: Q 1857 integer :: i 1858 1859 Q=0.D0 1860 IF (N>SIZE(X)) THEN 1861 call stop_application("QUANTILE1 [m_qtlmap_analyse] N ["//trim(str(N))//& 1862 '] to high [max='//trim(str(SIZE(X)))//']') 1863 END IF 1864 6010 DO 6011 I=1,N 1865 Q=Q+W(n,I,p)*X(I) 1866 6011 Continue 1867 RETURN 1868 END FUNCTION QUANTILE1
W
[ Top ] [ m_qtlmap_analyse ] [ Subroutines ]
NAME
W
DESCRIPTION
CALCUL DE W
INPUTS
n : i : p :
NOTES
use : MATH_QTLMAP_LOWERTAIL_BETA, MATH_QTLMAP_LOWERTAIL_BETA
SOURCE
1883 FUNCTION W(n,i,p) result(Wr) 1884 integer , intent(in) :: n 1885 integer , intent(in) :: i 1886 real (kind=dp) , intent(in) :: p 1887 real (kind=dp) :: Wr,TOL,X1,X2,associated,b 1888 real (kind=dp) :: p1,p2 1889 integer :: ifail 1890 1891 TOL=1.D-5 1892 X1=dble(i)/dble(n) 1893 X2=dble(i-1)/dble(n) 1894 associated=p*dble(n+1) 1895 b=(1.D0-p)*dble(n+1) 1896 ifail=1 1897 !CALL MATH_QTLMAP_G01EEF(X1,a,b,tol,p1,q1,pdf1,ifail) 1898 call MATH_QTLMAP_LOWERTAIL_BETA(X1, associated, b, TOL, p1, IFAIL) 1899 ifail=1 1900 !CALL MATH_QTLMAP_G01EEF(X2,a,b,tol,p2,q2,pdf2,ifail) 1901 call MATH_QTLMAP_LOWERTAIL_BETA(X2,associated,b,tol,p2,ifail) 1902 Wr=p1-p2 1903 RETURN 1904 END FUNCTION W