TOOLS
NAME
TOOLS
DESCRIPTION
Package tools : contains tools modules : mathematical and optimlization subroutines, lexer and parser subroutines 25/09/2009
m_qtlmap_tools
NAME
m_qtlmap_tools
DESCRIPTION
misc. routines.
NOTES
SEE ALSO
check_allocate
[ Top ] [ m_qtlmap_tools ] [ Subroutines ]
NAME
check_allocate
DESCRIPTION
check a iostat error code. print an error message and stop the programe if the iostat error /= 0
INPUTS
iostat : the iostat error code error_message : an error message
OUTPUTS
is_ok : boolean
RETURN
a double precision
SOURCE
377 subroutine check_allocate(iostat,error_message) 378 !ios_id where read the parameter 379 integer,intent(in)::iostat 380 !message to print if iostat /= 0 381 character(len=*),optional,intent(in)::error_message 382 383 if (iostat /= 0) then 384 print *, "Allocation error:[iostat:",iostat,"]" 385 if ( present (error_message) ) then 386 print *,'**'//error_message//'**' 387 endif 388 stop 'allocation error' 389 end if 390 end subroutine check_allocate
clean_lrt_solution
[ Top ] [ m_qtlmap_tools ] [ Subroutines ]
NAME
clean_lrt_solution
DESCRIPTION
release all sub data of TYPE_LRT_SOLUTION
INPUTS/OUTPUTS
lrtsol : the variable to release
SOURCE
1457 subroutine clean_lrt_solution(lrtsol) 1458 type(TYPE_LRT_SOLUTION) , intent(inout) :: lrtsol 1459 integer :: ip,jm 1460 1461 if ( lrtsol%nqtl < 0 ) return 1462 1463 if ( lrtsol%nqtl == 1 ) then 1464 if ( associated (lrtsol%lrt1) ) deallocate (lrtsol%lrt1) 1465 if ( associated (lrtsol%pater_eff) ) deallocate (lrtsol%pater_eff) 1466 if ( associated (lrtsol%mater_eff) ) deallocate (lrtsol%mater_eff) 1467 if ( associated (lrtsol%xlrp) ) deallocate (lrtsol%xlrp) 1468 if ( associated (lrtsol%xlrm) ) deallocate (lrtsol%xlrm) 1469 end if 1470 1471 if ( lrtsol%nqtl == 2 ) then 1472 if ( associated (lrtsol%lrt0_2) ) deallocate (lrtsol%lrt0_2) 1473 if ( associated (lrtsol%lrt1_2) ) deallocate (lrtsol%lrt1_2) 1474 if ( associated (lrtsol%pater_eff2) ) deallocate (lrtsol%pater_eff2) 1475 if ( associated (lrtsol%mater_eff2) ) deallocate (lrtsol%mater_eff2) 1476 if ( associated (lrtsol%xlrp2) ) deallocate (lrtsol%xlrp2) 1477 if ( associated (lrtsol%xlrm2) ) deallocate (lrtsol%xlrm2) 1478 end if 1479 1480 end subroutine clean_lrt_solution
copy
[ Top ] [ m_qtlmap_tools ] [ Subroutines ]
NAME
copy
DESCRIPTION
generic interface to copy a variable of qtlmap type
SOURCE
111 interface copy 112 module procedure copy_description_type 113 end interface copy
copy_description_type
[ Top ] [ m_qtlmap_tools ] [ Subroutines ]
NAME
copy_description_type
DESCRIPTION
copy the variable in_desc of type TYPE_INCIDENCE_SOLUTION to out_desc
INPUTS
in_desc : the variable to copy out_desc : the copy
SOURCE
1519 recursive subroutine copy_description_type(in_desc,out_desc) 1520 type(DESC_EFFECT_TYPE) , intent(in) :: in_desc 1521 type(DESC_EFFECT_TYPE) , intent(out) :: out_desc 1522 integer :: i 1523 1524 out_desc%name = in_desc%name 1525 out_desc%isVar = in_desc%isVar 1526 out_desc%start = in_desc%start 1527 out_desc%end = in_desc%end 1528 out_desc%haveSubDesc = in_desc%haveSubDesc 1529 1530 if ( in_desc%haveSubDesc ) then 1531 allocate ( out_desc%listSubDesc(size(in_desc%listSubDesc)) ) 1532 do i=1,size(in_desc%listSubDesc) 1533 call copy_description_type(in_desc%listSubDesc(i),out_desc%listSubDesc(i)) 1534 end do 1535 end if 1536 1537 end subroutine copy_description_type
dget_str
[ Top ] [ m_qtlmap_tools ] [ Subroutines ]
NAME
dget_str
DESCRIPTION
conversion of a double precision into a character word
INPUTS
value : the word
RETURN
a double precision
SOURCE
563 function dget_str(value) result(res) 564 real (kind=dp) ,intent(in) :: value 565 character(len=20) :: res 566 character(len=LEN_DEF) :: str_value,c,sign_value,tempvs,tempvs2,temppoint 567 integer :: ascii_zero,reste,decimal 568 integer :: val_int,io 569 real(kind=dp) :: v 570 !real :: v 571 logical :: have_dec 572 573 res='' 574 write (res,fmt='(f13.2)',iostat=io) value 575 if ( io /= 0 ) then 576 res='<NO-DEF>' 577 else 578 res=adjustl(res) 579 end if 580 return 581 582 end function dget_str
EXTRACT
[ Top ] [ m_qtlmap_tools ] [ Subroutines ]
NAME
EXTRACT
DESCRIPTION
Get a subset of a character line
INPUTS
string : the line start : start index of the subset finish : end index of the subset
RETURN
the subset of the line
SOURCE
696 function EXTRACT (string,start,finish) result(out) 697 character(len=*) , intent(in) :: string 698 integer , intent(in) ,optional :: start 699 integer , intent(in) ,optional :: finish 700 701 character(len=len(string)) :: out 702 integer :: s,e 703 704 s=1 705 e=len(string) 706 707 if (present(start)) then 708 if (start > 1 .and. start <= len(string)) then 709 s = start 710 end if 711 end if 712 713 if (present(finish)) then 714 if (finish >= start .and. finish <= len(string)) then 715 e = finish 716 end if 717 end if 718 719 out=trim(string(s:e)) 720 721 end function
fget_str
[ Top ] [ m_qtlmap_tools ] [ Subroutines ]
NAME
fget_str
DESCRIPTION
conversion of a a float into a character word
INPUTS
value : the word
RETURN
a float
SOURCE
540 function fget_str(value) result(res) 541 real ,intent(in) :: value 542 character(len=20) :: res 543 real (kind=dp) :: v 544 545 v = value 546 res = dget_str(v) 547 548 end function fget_str
file_exist
[ Top ] [ m_qtlmap_tools ] [ Subroutines ]
NAME
file_exist
DESCRIPTION
Test if a local file exist (the programe stop if not)
INPUTS
file_in : the pathname file
SOURCE
175 subroutine file_exist(file_in) 176 character(len=*),intent(inout) ::file_in 177 logical::existe 178 inquire( file=file_in, exist=existe) 179 180 if(.not. existe) then 181 call stop_application('The file ['//trim(file_in)//'] does not exist.') 182 endif 183 end subroutine file_exist
GET
[ Top ] [ m_qtlmap_tools ] [ Subroutines ]
NAME
GET
DESCRIPTION
Read a line from an unit file
INPUTS
unit : the unit record maxlen : maximum length of the line
INPUTS/OUTPUTS
string : the line read
OUTPUTS
iostat : iostat error
SOURCE
739 subroutine GET (unit,string,maxlen,iostat) 740 integer , intent(in) :: unit 741 character(len=LEN_LINE) , intent(inout) :: string 742 integer , intent(in),optional :: maxlen 743 integer , intent(out),optional :: iostat 744 745 character(len=LEN_LINE) :: word 746 integer :: eof 747 748 eof = 0 749 string='' 750 read (unit,iostat=eof,fmt="(a)/") string 751 if ( present(iostat) ) iostat = eof 752 753 end subroutine GET
get_dx
[ Top ] [ m_qtlmap_tools ] [ Subroutines ]
NAME
get_dx
DESCRIPTION
Get the position in Morgan
INPUTS
index_chr : index of chromosome ix : the position index
RETURN
the position in Morgan
NOTES
SOURCE
983 function get_dx(index_chr,ix) result(dx) 984 integer, intent(in) :: index_chr,ix 985 real(kind=dp) :: dx 986 987 if ( index_chr > nchr ) then 988 print *,'ind:',index_chr 989 print *,'dev error: get_ilong can not access to index chr' 990 stop 991 end if 992 993 dx=posi(index_chr,1)+(dble(ix)/BASE_STEP) 994 ! print *,"ix:",ix,"dx:",dx 995 return 996 997 end function get_dx
get_ilong
[ Top ] [ m_qtlmap_tools ] [ Subroutines ]
NAME
get_ilong
DESCRIPTION
Get the number of sampled position among the selected chromosome
INPUTS
index_chr : index of chromosome RESULT The chromosome size in number of position sampled
NOTES
SEE ALSO
m_qtlmap_data/nmk m_qtlmap_data/posi
SOURCE
950 function get_ilong(index_chr) result(ilong) 951 integer, intent(in) :: index_chr 952 integer :: ilong 953 954 if ( index_chr > nchr ) then 955 print *,'ind:',index_chr 956 print *,'dev error: get_ilong can not access to index chr' 957 stop 958 end if 959 ilong=nint(BASE_STEP*(posi(index_chr,nmk(index_chr))-posi(index_chr,1))) 960 return 961 962 end function get_ilong
get_ind_pheno
[ Top ] [ m_qtlmap_tools ] [ Subroutines ]
NAME
get_ind_pheno
DESCRIPTION
Get the corresponding internal identifiant number of a coded initial allele string
INPUTS
* c : the chromosome * value : Value of the string allele RESULT The identifiant value of the value allele
EXAMPLE
NOTES
SEE ALSO
m_qtlmap_data/value_pheno m_qtlmap_data/VAL_MIN_INDEX_PHENO
SOURCE
1237 function get_ind_pheno(c,value) result(res) 1238 integer ,intent(in) :: c 1239 character(len=LEN_DEF) ,intent(in) :: value 1240 integer(kind=KIND_PHENO) :: res 1241 integer :: i 1242 1243 do i=1,nb_value_pheno(c) 1244 if ( value_pheno(c,i) == value ) then 1245 res = VAL_MIN_INDEX_PHENO + i + 1 1246 return 1247 end if 1248 end do 1249 write (0,*) "Devel error : get_ind_pheno * none value correspond [",trim(value),"]" 1250 stop 1251 end function get_ind_pheno
get_int
[ Top ] [ m_qtlmap_tools ] [ Subroutines ]
NAME
get_int
DESCRIPTION
generic interface to convert a word into an integer
SOURCE
147 interface get_int 148 module procedure get_int_char 149 end interface get_int
get_int_char
[ Top ] [ m_qtlmap_tools ] [ Subroutines ]
NAME
get_int_char
DESCRIPTION
conversion of a character word in an integer using the fortran read routine
INPUTS
string_to_convert : the word
OUTPUTS
is_ok : boolean
RETURN
an integer
SOURCE
462 function get_int_char(string_to_convert,is_ok) result(int_value) 463 character(len=*),intent(in) :: string_to_convert 464 logical,optional,intent(inout) :: is_ok 465 integer :: int_value 466 character(len=LEN_BUFFER_WORD) :: word_token 467 integer :: ios 468 469 word_token = trim(string_to_convert)!CHAR(string_to_convert) 470 read (word_token,FMT_INT,IOSTAT=ios) int_value 471 472 if ( present(is_ok) ) then 473 is_ok = ( ios == 0 ) 474 endif 475 476 if ( ios /= 0 ) int_value = INT_NOT_DEFINED 477 word_token ='' 478 return 479 480 end function get_int_char
get_long_step_morgan
[ Top ] [ m_qtlmap_tools ] [ Subroutines ]
NAME
get_long_step_morgan
DESCRIPTION
Get the the step in Morgan
SOURCE
809 function get_long_step_morgan() result(l) 810 real(kind=dp) :: l 811 812 l = real(pas) / BASE_STEP 813 814 end function get_long_step_morgan
get_maxnbgenotypedam
[ Top ] [ m_qtlmap_tools ] [ Subroutines ]
NAME
get_maxnbgenotypedam
DESCRIPTION
get the maximum number of genotype (all dams confused) found by the haplotype on the genome wide RESULTS maxng : the variable to release
SOURCE
1601 function get_maxnbgenotypedam() result(maxng) 1602 integer :: maxng 1603 1604 maxng = maxval(ngenom(:,nmp(np+1)+1)) 1605 1606 end function get_maxnbgenotypedam
get_maxnpo
[ Top ] [ m_qtlmap_tools ] [ Subroutines ]
NAME
get_npo
DESCRIPTION
Getting the maximum npo RESULT maximum number of position on a chromosome
EXAMPLE
NOTES
SEE ALSO
m_qtlmap_data/get_ilong
SOURCE
1071 function get_maxnpo() result(maxnpo) 1072 integer :: chr,npo,maxnpo 1073 1074 maxnpo=0 1075 1076 do chr=1,nchr 1077 npo = get_npo(chr) 1078 maxnpo = max(maxnpo,npo) 1079 end do 1080 1081 return 1082 1083 end function get_maxnpo
get_npo
[ Top ] [ m_qtlmap_tools ] [ Subroutines ]
NAME
get_npo
DESCRIPTION
Getting the number of postion depending the chromosome , step analysis and ilong (size of the chromosome)
INPUTS
* index_chr -- Chromosome number RESULT The number of position to computing Probabilities of transmission and likelihood
EXAMPLE
NOTES
SEE ALSO
DATA/m_qtlmap_data/get_ilong
SOURCE
1046 function get_npo(index_chr) result(npo) 1047 integer, intent(in) :: index_chr 1048 integer :: npo 1049 npo = (get_ilong(index_chr) / pas)+1 1050 return 1051 1052 end function get_npo
get_pheno
[ Top ] [ m_qtlmap_tools ] [ Subroutines ]
NAME
get_pheno
DESCRIPTION
Get the corresponding identifiant number of a coded initial allele string
INPUTS
* value -- Value of the allele string RESULT The string value of the allele
EXAMPLE
NOTES
SEE ALSO
m_qtlmap_data/value_pheno m_qtlmap_data/VAL_MIN_INDEX_PHENO
SOURCE
1201 function get_pheno(c,value) result(res) 1202 integer ,intent(in) :: c 1203 integer(kind=KIND_PHENO) ,intent(in) :: value 1204 character(len=LEN_DEF) :: res 1205 integer :: ind,i 1206 1207 ind = value - VAL_MIN_INDEX_PHENO + 1 1208 if ( ind < 1 .or. ind > nb_value_pheno(c)) then 1209 write (0,*) "Devel error : get_pheno * none index correspond [",value,"]" 1210 print *, (trim(value_pheno(c,i))//" ",i=1,nb_value_pheno(c)) 1211 stop 1212 end if 1213 res = value_pheno(c,ind) 1214 end function get_pheno
get_pos
[ Top ] [ m_qtlmap_tools ] [ Subroutines ]
NAME
get_pos
DESCRIPTION
Get the index of the sampled position
INPUTS
posInMorgan : position in Morgan
RETURN
the index of sampled position
NOTES
SOURCE
1016 function get_pos(posInMorgan) result(npo) 1017 real(kind=dp), intent(in) :: posInMorgan 1018 integer :: npo 1019 1020 npo=int(BASE_STEP*posInMorgan) 1021 npo=nint(real(npo)/real(PAS)) 1022 1023 return 1024 1025 end function get_pos
get_real
[ Top ] [ m_qtlmap_tools ] [ Subroutines ]
NAME
get_real
DESCRIPTION
generic interface to convert a word into a float
SOURCE
136 interface get_real 137 module procedure get_real_char 138 end interface get_real
get_real_char
[ Top ] [ m_qtlmap_tools ] [ Subroutines ]
NAME
get_real_char
DESCRIPTION
conversion of a character word in a double precision using the fortran read routine
INPUTS
string_to_convert : the word
OUTPUTS
is_ok : boolean
RETURN
a double precision
SOURCE
407 function get_real_char(string_to_convert,is_ok) result(real_value) 408 character(len=*),intent(in) :: string_to_convert 409 logical,optional,intent(inout) :: is_ok 410 real (kind=dp) :: real_value 411 character(len=LEN_BUFFER_WORD) :: word_token 412 integer :: ios 413 character(len=LEN_DEF) :: c,string 414 logical :: is_int 415 416 string = trim(string_to_convert) 417 is_int=.true. 418 419 do while ( string /='' ) 420 c = EXTRACT (string , 1, 1) 421 if (c == '.') then 422 is_int=.false. 423 exit 424 end if 425 string = REMOVE(string,1,1) 426 end do 427 428 if ( .not. is_int ) then 429 word_token = string_to_convert!CHAR(string_to_convert) 430 read (word_token,FMT_REAL,IOSTAT=ios) real_value 431 if ( present(is_ok) ) then 432 is_ok = ( ios == 0 ) 433 end if 434 if ( ios /= 0 ) real_value = REAL_NOT_DEFINED 435 else 436 if (present(is_ok)) then 437 real_value = get_int(string_to_convert,is_ok) 438 else 439 real_value = get_int(string_to_convert) 440 end if 441 end if 442 string='' 443 return 444 445 end function get_real_char
iget_str
[ Top ] [ m_qtlmap_tools ] [ Subroutines ]
NAME
iget_str
DESCRIPTION
conversion of an integer into a character word
INPUTS
value_int : the word
RETURN
an integer
SOURCE
495 function iget_str(value_int) result(res) 496 integer,intent(in) :: value_int 497 character(len=20) :: res 498 character(len=LEN_DEF) :: str_value,sign_value,tempvs,tempvs2 499 integer :: v,reste,ascii_zero 500 501 if (value_int<0) then 502 sign_value='-' 503 v = -value_int 504 else 505 sign_value='' 506 v = value_int 507 end if 508 509 ascii_zero = iachar('0') 510 str_value='' 511 do while (v /=0 ) 512 reste = mod(v,10) 513 tempvs=achar(ascii_zero+reste) 514 tempvs2=trim(str_value) 515 str_value = trim(tempvs)//trim(tempvs2) 516 v = v / 10 517 end do 518 519 str_value = trim(sign_value)//trim(str_value) 520 if (str_value=='' .or. str_value=='-') str_value='0' 521 res = trim(str_value) 522 523 return 524 525 end function iget_str
init_pheno
[ Top ] [ m_qtlmap_tools ] [ Subroutines ]
NAME
init_pheno
DESCRIPTION
allocate arrays to manage the internal encoded format of allele : value_pheno,nb_value_pheno
INPUTS
nch : number of chromosome to manage
NOTES
SEE ALSO
SOURCE
1099 subroutine init_pheno(nch) 1100 integer ,intent(in) :: nch 1101 allocate (value_pheno(nch,256)) 1102 allocate (nb_value_pheno(nch)) 1103 value_pheno="" 1104 nb_value_pheno=0 1105 1106 end subroutine init_pheno
new
[ Top ] [ m_qtlmap_tools ] [ Subroutines ]
NAME
new
DESCRIPTION
generic interface to allocate internal array of a qtlmap variable with a complex type
SOURCE
87 interface new 88 module procedure new_lrt_solution 89 end interface new
new_lrt_solution
[ Top ] [ m_qtlmap_tools ] [ Subroutines ]
NAME
new_lrt_solution
DESCRIPTION
initialize a variable of type TYPE_LRT_SOLUTION
INPUTS/OUTPUTS
lrtsol : the variable to initialize
SOURCE
1363 subroutine new_lrt_solution(nqtl,lrtsol) 1364 integer , intent(in) :: nqtl 1365 type(TYPE_LRT_SOLUTION) , intent(inout) :: lrtsol 1366 integer :: npo,i,iq,ip,jm 1367 ! solution under hypothesis nqtl 1368 lrtsol%nqtl=nqtl 1369 1370 ! begin index <--0 1371 allocate (lrtsol%lrtmax(0:nqtl-1)) 1372 allocate (lrtsol%nxmax(0:nqtl-1)) 1373 allocate (lrtsol%chrmax(0:nqtl-1)) 1374 1375 lrtsol%lrtmax = 0.d0 1376 lrtsol%nxmax = 0 1377 lrtsol%chrmax = 0 1378 1379 if ( nqtl == 1 .or. nqtl == 2 ) then 1380 npo = get_maxnpo() 1381 if ( nqtl == 1) then 1382 allocate (lrtsol%lrt1(nchr,npo)) 1383 allocate (lrtsol%pater_eff(nchr,np,npo)) 1384 allocate (lrtsol%mater_eff(nchr,nm,npo)) 1385 allocate (lrtsol%xlrp(nchr,np,npo)) 1386 allocate (lrtsol%xlrm(nchr,nm,npo)) 1387 lrtsol%pater_eff=0.d0 1388 lrtsol%mater_eff=0.d0 1389 lrtsol%xlrp=0.d0 1390 lrtsol%xlrm=0.d0 1391 end if 1392 1393 if ( nqtl == 2 ) then 1394 allocate (lrtsol%lrt0_2(nchr,nchr,npo,npo)) 1395 allocate (lrtsol%lrt1_2(nchr,nchr,npo,npo)) 1396 allocate (lrtsol%pater_eff2(nchr,nchr,np,npo,npo,2)) 1397 allocate (lrtsol%mater_eff2(nchr,nchr,nm,npo,npo,2)) 1398 allocate (lrtsol%xlrp2(nchr,nchr,np,npo,npo)) 1399 allocate (lrtsol%xlrm2(nchr,nchr,nm,npo,npo)) 1400 lrtsol%lrt0_2=0.d0 1401 lrtsol%lrt1_2=0.d0 1402 lrtsol%pater_eff2=0.d0 1403 lrtsol%mater_eff2=0.d0 1404 lrtsol%xlrp2=0.d0 1405 lrtsol%xlrm2=0.d0 1406 end if 1407 end if 1408 end subroutine new_lrt_solution
next_word
[ Top ] [ m_qtlmap_tools ] [ Subroutines ]
NAME
next_word
DESCRIPTION
generic interface to read a word from line (parsing context)
SOURCE
158 interface next_word 159 module procedure next_word_char 160 end interface next_word
next_word_char
[ Top ] [ m_qtlmap_tools ] [ Subroutines ]
NAME
next_word_char
DESCRIPTION
get the last word at the right side of a line and remove it. INPUTS/OUTPUT line_to_parse : the line.
OUTPUTS
is_ok : a flag error
RETURN
the word
SOURCE
308 FUNCTION next_word_char(line_to_parse,is_ok) result(w_token) 309 310 character(len=*) ,intent(inout) :: line_to_parse 311 logical, intent(inout),optional :: is_ok 312 character(len=LEN_BUFFER_WORD) :: w_token 313 integer :: nd 314 315 nd = 1 316 w_token = '' 317 line_to_parse=adjustl(line_to_parse) 318 319 if ( trim(line_to_parse)=='' ) then 320 if ( present(is_ok) ) is_ok = .false. 321 return 322 end if 323 324 325 do while (nd == 1 .and. w_token == '') 326 !! char(9) pour les tabulation..... 327 nd=SCAN(line_to_parse," ;"//char(9)) 328 329 if ( nd == 0) then 330 w_token = trim(line_to_parse) 331 else 332 w_token = trim(EXTRACT(line_to_parse,1,nd-1)) 333 end if 334 335 336 if (w_token == '') then 337 line_to_parse = REMOVE(line_to_parse,1,nd) 338 end if 339 end do 340 341 if (nd == 0) then ! END OF LINE 342 w_token = trim(EXTRACT(line_to_parse,1)) 343 line_to_parse = '' 344 if ( present(is_ok) ) is_ok = .true. 345 return 346 end if 347 348 if (nd > 0) then ! at least one non-terminator character in the word 349 w_token = trim(EXTRACT(line_to_parse,1,nd-1)) 350 if ( present(is_ok) ) is_ok = .true. 351 else 352 w_token = '' 353 if ( present(is_ok) ) is_ok = .false. 354 end if 355 356 357 line_to_parse = adjustl(REMOVE(line_to_parse,1,nd)) 358 return 359 end function next_word_char
parse_real_array
[ Top ] [ m_qtlmap_tools ] [ Subroutines ]
NAME
parse_real_array
DESCRIPTION
get a set of real value from a line
INPUTS
nvalue : number of value INPUTS/OUTPUT line_to_parse : the line
OUTPUTS
tabV : a vector of double precision is_ok : a flag error
SOURCE
246 subroutine parse_real_array(line_to_parse,nvalue,tabV,is_ok) 247 character(len=*),intent(inout) :: line_to_parse 248 integer, intent(in) :: nvalue 249 real(kind=dp) ,dimension(nvalue),intent(out) :: tabV 250 logical, intent(inout),optional :: is_ok 251 character(len=LEN_BUFFER_WORD) :: w_token 252 integer :: i,ns,ne 253 character(len=1000) :: endline 254 255 if ( trim(line_to_parse)=='' ) then 256 if ( present(is_ok) ) is_ok = .false. 257 return 258 end if 259 260 ns=index(line_to_parse,"[") 261 if ( ns == 0 ) then 262 if ( present(is_ok) ) is_ok = .false. 263 return 264 end if 265 ne=index(line_to_parse,"]") 266 267 if ( ne == 0 ) then 268 if ( present(is_ok) ) is_ok = .false. 269 return 270 end if 271 272 endline = line_to_parse(ne+1:) 273 line_to_parse = line_to_parse(ns+1:ne-1) 274 275 do i=1,nvalue 276 w_token=trim(next_word_char(line_to_parse,is_ok)) 277 if (is_ok) tabV(i)=get_real(w_token,is_ok) 278 if ( .not. is_ok ) then 279 if ( present(is_ok) ) is_ok = .false. 280 return 281 end if 282 end do 283 284 if ( trim(line_to_parse) /= '' ) then 285 if ( present(is_ok) ) is_ok = .false. 286 return 287 end if 288 289 line_to_parse = endline 290 291 end subroutine parse_real_array
release
[ Top ] [ m_qtlmap_tools ] [ Subroutines ]
NAME
release
DESCRIPTION
generic interface to release internal array of a qtlmap variable with a complex type
SOURCE
98 interface release 99 module procedure release_lrt_solution, & 100 release_incidence_solution, & 101 release_description_type 102 end interface release
release_data
[ Top ] [ m_qtlmap_tools ] [ Subroutines ]
NAME
release_data
DESCRIPTION
release all arrays from the module m_qtlmap_data
SOURCE
1260 subroutine release_data 1261 integer :: ic 1262 1263 deallocate (repro) 1264 deallocate (gpere) 1265 deallocate (gmere) 1266 deallocate (pere) 1267 deallocate (mere) 1268 deallocate (ngmgp) 1269 deallocate (nrgm) 1270 deallocate (reppere) 1271 deallocate (repmere) 1272 deallocate (animal) 1273 deallocate (nmp) 1274 deallocate (ndm) 1275 deallocate (femelle) 1276 deallocate (repfem) 1277 1278 deallocate (nmk) 1279 deallocate (mark) 1280 if (allocated(ch)) deallocate (ch) 1281 deallocate (posi) 1282 deallocate (posim) 1283 deallocate (posif) 1284 deallocate (rm) 1285 deallocate (rf) 1286 deallocate (absi) 1287 1288 deallocate (numero) 1289 deallocate (pheno) 1290 deallocate (value_pheno) 1291 deallocate (nb_value_pheno) 1292 deallocate (corregp) 1293 deallocate (corregm) 1294 deallocate (correr) 1295 deallocate (correp) 1296 deallocate (correm) 1297 deallocate (corred) 1298 deallocate (nall) 1299 deallocate (presentg) 1300 deallocate (alleles) 1301 deallocate (pc_all) 1302 1303 if (allocated(bete)) deallocate (bete) 1304 deallocate (modele) 1305 deallocate (carac) 1306 deallocate (natureY) 1307 deallocate (y) 1308 deallocate (cd) 1309 deallocate (xmut) 1310 deallocate (sigt) 1311 deallocate (ydiscretord) 1312 if (allocated (ycategorial)) deallocate (ycategorial) 1313 deallocate (indicemod) 1314 deallocate (nmod) 1315 deallocate (seuil) 1316 deallocate (prop) 1317 deallocate (presentc) 1318 deallocate (ndelta) 1319 deallocate (corperf) 1320 deallocate (namefix) 1321 deallocate (nlev) 1322 deallocate (listelev) 1323 deallocate (nivx) 1324 deallocate (namecov) 1325 deallocate (covar) 1326 deallocate (chromo) 1327 1328 deallocate (estime) 1329 deallocate (nmumest) 1330 deallocate (namest) 1331 if (allocated(iam)) deallocate (iam) 1332 deallocate (estfem) 1333 1334 if (allocated(h2)) deallocate (h2) 1335 if (allocated(RhoP)) deallocate (RhoP) 1336 if (allocated(RhoG)) deallocate (RhoG) 1337 if (allocated(ue)) deallocate (ue) 1338 if (allocated(soglia)) deallocate (soglia) 1339 if (allocated(freqmodsimultrait)) deallocate (freqmodsimultrait) 1340 if (allocated(nbmodsimultrait)) deallocate (nbmodsimultrait) 1341 if (allocated(posiqtl)) deallocate (posiqtl) 1342 if (allocated(chrqtl)) deallocate (chrqtl) 1343 if (allocated(nbmodsimultrait)) deallocate (xlim) 1344 1345 if (associated (listModelTrait)) then 1346 do ic=1,ncar 1347 call release_model_trait(listModelTrait(ic)) 1348 end do 1349 deallocate (listModelTrait) 1350 end if 1351 1352 end subroutine release_data
release_description_type
[ Top ] [ m_qtlmap_tools ] [ Subroutines ]
NAME
release_description_type
DESCRIPTION
release a variable of type DESC_EFFECT_TYPE
INPUTS/OUTPUTS
inout_desc : the variable to release
SOURCE
1548 recursive subroutine release_description_type(inout_desc) 1549 type(DESC_EFFECT_TYPE) , intent(inout) :: inout_desc 1550 integer :: i 1551 1552 if ( inout_desc%haveSubDesc ) then 1553 do i=1,size(inout_desc%listSubDesc) 1554 call release_description_type(inout_desc%listSubDesc(i)) 1555 end do 1556 deallocate (inout_desc%listSubDesc) 1557 end if 1558 1559 end subroutine release_description_type
release_incidence_solution
[ Top ] [ m_qtlmap_tools ] [ Subroutines ]
NAME
release_incidence_solution
DESCRIPTION
release a variable of type TYPE_INCIDENCE_SOLUTION
INPUTS/OUTPUTS
obj : the variable to release
SOURCE
1491 subroutine release_incidence_solution(obj) 1492 type(TYPE_INCIDENCE_SOLUTION) , intent(inout) :: obj 1493 1494 if ( associated (obj%sig) ) deallocate (obj%sig) 1495 if ( associated (obj%unknown_dam_sig) ) deallocate (obj%unknown_dam_sig) 1496 if ( associated (obj%paramaterValue) ) deallocate (obj%paramaterValue) 1497 if ( associated (obj%parameterVecsol) ) deallocate (obj%parameterVecsol) 1498 if ( associated (obj%eqtl_print) ) deallocate (obj%eqtl_print) 1499 if ( associated (obj%parameterPrecis) ) deallocate (obj%parameterPrecis) 1500 if ( associated(obj%rhoi) ) deallocate (obj%rhoi) 1501 1502 if ( associated (obj%groupeName) ) deallocate (obj%groupeName) 1503 if ( associated (obj%nbParameterGroup) ) deallocate (obj%nbParameterGroup) 1504 if ( associated (obj%parameterName) ) deallocate (obj%parameterName) 1505 if (associated(obj%qtl_groupeName)) deallocate (obj%qtl_groupeName) 1506 1507 end subroutine release_incidence_solution
release_lrt_solution
[ Top ] [ m_qtlmap_tools ] [ Subroutines ]
NAME
release_lrt_solution
DESCRIPTION
release a variable of type TYPE_LRT_SOLUTION
INPUTS/OUTPUTS
lrtsol : the variable to release
SOURCE
1419 subroutine release_lrt_solution(lrtsol) 1420 type(TYPE_LRT_SOLUTION) , intent(inout) :: lrtsol 1421 integer :: ip,jm 1422 1423 if ( lrtsol%nqtl < 0 ) return 1424 1425 if ( associated (lrtsol%lrtmax) ) deallocate (lrtsol%lrtmax) 1426 if ( associated (lrtsol%chrmax) ) deallocate (lrtsol%chrmax) 1427 if ( associated (lrtsol%nxmax) ) deallocate (lrtsol%nxmax) 1428 1429 if ( lrtsol%nqtl == 1 ) then 1430 if ( associated (lrtsol%lrt1) ) deallocate (lrtsol%lrt1) 1431 if ( associated (lrtsol%pater_eff) ) deallocate (lrtsol%pater_eff) 1432 if ( associated (lrtsol%mater_eff) ) deallocate (lrtsol%mater_eff) 1433 if ( associated (lrtsol%xlrp) ) deallocate (lrtsol%xlrp) 1434 if ( associated (lrtsol%xlrm) ) deallocate (lrtsol%xlrm) 1435 end if 1436 1437 if ( lrtsol%nqtl == 2 ) then 1438 if ( associated (lrtsol%lrt0_2) ) deallocate (lrtsol%lrt0_2) 1439 if ( associated (lrtsol%lrt1_2) ) deallocate (lrtsol%lrt1_2) 1440 if ( associated (lrtsol%pater_eff2) ) deallocate (lrtsol%pater_eff2) 1441 if ( associated (lrtsol%mater_eff2) ) deallocate (lrtsol%mater_eff2) 1442 if ( associated (lrtsol%xlrp2) ) deallocate (lrtsol%xlrp2) 1443 if ( associated (lrtsol%xlrm2) ) deallocate (lrtsol%xlrm2) 1444 end if 1445 1446 end subroutine release_lrt_solution
release_model_trait
[ Top ] [ m_qtlmap_tools ] [ Subroutines ]
NAME
release_model_trait
DESCRIPTION
release a variable of type model_trait
INPUTS/OUTPUTS
inout_desc : the variable to release
SOURCE
1570 recursive subroutine release_model_trait(inout_desc) 1571 type(model_trait) , intent(inout) :: inout_desc 1572 1573 if ( associated(inout_desc%nbint) ) then 1574 deallocate (inout_desc%nbint) 1575 end if 1576 1577 if ( associated(inout_desc%indexFixedEffect) ) then 1578 deallocate (inout_desc%indexFixedEffect) 1579 end if 1580 1581 if ( associated(inout_desc%indexCovariate) ) then 1582 deallocate (inout_desc%indexCovariate) 1583 end if 1584 1585 if ( associated(inout_desc%indexFixedEffectWithInteraction) ) then 1586 deallocate (inout_desc%indexFixedEffectWithInteraction) 1587 end if 1588 1589 end subroutine release_model_trait
REMOVE
[ Top ] [ m_qtlmap_tools ] [ Subroutines ]
NAME
REMOVE
DESCRIPTION
Remove a subset of a line
INPUTS
string : the line start : start index of the subset finish : end index of the subset
RETURN
the new line without the subset
SOURCE
646 function REMOVE (string,start,finish) result(out) 647 character(len=*) , intent(in) :: string 648 integer , intent(in) ,optional :: start 649 integer , intent(in) ,optional :: finish 650 651 character(len=len(string)) :: out 652 integer :: s,e 653 654 s = 1 655 e = len(string) 656 if (present(start)) then 657 if (start > 1 .and. start <= len(string)) then 658 s = start 659 end if 660 end if 661 662 if (present(finish)) then 663 if (finish >= start .and. finish <= len(string)) then 664 e = finish 665 end if 666 end if 667 668 out='' 669 670 if ( s > 1) then 671 out=trim(string(1:s)) 672 end if 673 674 if ( e < len(string) ) then 675 out = trim(out) // string(e+1:len(string)) 676 end if 677 678 end function REMOVE
set_absi
[ Top ] [ m_qtlmap_tools ] [ Subroutines ]
NAME
set_absi
DESCRIPTION
build the global array absi. absi is a corresponding vector with a sampled positon. absi(ch,iposi)=<position in morgan of iposi>, iposi is an integer 1<=iposi<=get_ilong(ch)
SOURCE
826 subroutine set_absi() 827 integer :: i,chr,l,n,ix,ilk,nlong 828 real :: v,lastv,lastvm,lastvf 829 logical :: vOk,vmOk,vfOk 830 831 integer ,allocatable , dimension(:,:) :: corr_average_absi,corr_mal_absi,corr_femal_absi 832 833 allocate (corr_average_absi(nchr,maxval(nmk))) 834 allocate (corr_mal_absi(nchr,maxval(nmk))) 835 allocate (corr_femal_absi(nchr,maxval(nmk))) 836 corr_average_absi=0 837 corr_mal_absi=0 838 corr_femal_absi=0 839 840 if ( nchr == 0 ) then 841 call stop_application("Devel error: Can not set absci before the map structure!") 842 end if 843 844 allocate (absi(nchr,get_maxnpo())) 845 846 do chr=1,nchr 847 l = get_ilong(chr) 848 n=0 849 do ix=0,l,pas 850 n=n+1 851 absi(chr,n)=get_dx(chr,ix) 852 end do 853 854 nlong=n 855 do ilk=1,nmk(chr) 856 lastv = 10000.d0;lastvm = 10000.d0;lastvf = 10000.d0 857 vOk = .false. 858 vmOk = .false. 859 vfOk = .false. 860 do n=1,nlong 861 !! average map 862 863 v = ( posi(chr,ilk) - absi(chr,n) ) / posi(chr,ilk) 864 865 if ( v <= 0.0 .or. n == nlong ) then 866 867 if ( v < 0.d0 ) then ! maybe the last value was nearest.... 868 if ( abs(v) < abs(lastv) ) then 869 corr_average_absi(chr,ilk)=n 870 else 871 corr_average_absi(chr,ilk)=n-1 872 end if 873 else 874 corr_average_absi(chr,ilk)=n 875 end if 876 vOk = .true. 877 end if 878 lastv = v 879 !! mal map 880 881 v = ( posim(chr,ilk) - absi(chr,n) ) / posim(chr,ilk) 882 883 if ( v <= 0.0 .or. n == nlong ) then 884 885 if ( v < 0.d0 ) then ! maybe the last value was nearest.... 886 if ( abs(v) < abs(lastvm) ) then 887 corr_mal_absi(chr,ilk)=n 888 else 889 corr_mal_absi(chr,ilk)=n-1 890 end if 891 else 892 corr_mal_absi(chr,ilk)=n 893 end if 894 vmOk = .true. 895 end if 896 lastvm = v 897 898 !! femal map 899 900 v = ( posif(chr,ilk) - absi(chr,n) ) / posif(chr,ilk) 901 902 if ( v <= 0.0 .or. n == nlong ) then 903 904 if ( v < 0.d0 ) then ! maybe the last value was nearest.... 905 if ( abs(v) < abs(lastvf) ) then 906 corr_femal_absi(chr,ilk)=n 907 else 908 corr_femal_absi(chr,ilk)=n-1 909 end if 910 else 911 corr_femal_absi(chr,ilk)=n 912 end if 913 vfOk = .true. 914 end if 915 lastvf = v 916 917 918 if ( vOk .and. vmOk .and. vfOk ) exit 919 920 end do 921 call log_mess("CH "//str(chr)//" m"//str(ilk)//":"//& 922 str(corr_average_absi(chr,ilk))//" "//& 923 str(corr_mal_absi(chr,ilk))//" "//& 924 str(corr_femal_absi(chr,ilk)),DEBUG_DEF); 925 end do 926 end do 927 928 deallocate (corr_average_absi) 929 deallocate (corr_mal_absi) 930 deallocate (corr_femal_absi) 931 ! stop 932 end subroutine set_absi
set_base_and_step
[ Top ] [ m_qtlmap_tools ] [ Subroutines ]
NAME
set_base_and_step
DESCRIPTION
Set the step (variable PAS) and the base dimension (variable BASE_STEP) according a Morgan step.
INPUTS
stepMorgan : the step in Morgan
SOURCE
766 subroutine set_base_and_step(stepMorgan) 767 character(len=*) ,intent(in) :: stepMorgan 768 integer :: i,s 769 character(len=LEN_W) :: buf 770 771 buf="" 772 buf=adjustl(stepMorgan) 773 s=len(trim(buf)) 774 !on cherche le point..... 775 do i=1,s 776 if (buf(i:i) == '.' ) then 777 exit 778 end if 779 end do 780 ! pas de point... 781 if ( i > s ) then 782 BASE_STEP = 1 783 PAS = get_int(buf) 784 else 785 if ( (s-i) >= 10 ) then 786 call stop_application("** Bad definition of step **") 787 end if 788 BASE_STEP = 10**(s-i) 789 buf=trim(buf(:i-1))//trim(buf(i+1:)) 790 do while ( buf(1:1)=='0' ) 791 buf=trim(buf(2:)) 792 end do 793 PAS = get_int(buf) 794 end if 795 796 call log_mess("BASE STEP : "//trim(str(BASE_STEP)),DEBUG_DEF) 797 call log_mess(" STEP : "//trim(str(PAS)),DEBUG_DEF) 798 799 end subroutine set_base_and_step
set_pheno
[ Top ] [ m_qtlmap_tools ] [ Subroutines ]
NAME
set_pheno
DESCRIPTION
Setting a new identifiant number allele with allele string value or get the value of the identifiant if this value exist. * Reallocation or Allocation (first call) of the value_pheno for each call
INPUTS
* value -- Value of the allele string RESULT The idenfiant allele number
EXAMPLE
NOTES
SEE ALSO
m_qtlmap_data/value_pheno m_qtlmap_data/VAL_MIN_INDEX_PHENO m_qtlmap_data/VAL_MAX_INDEX_PHENO
SOURCE
1131 function set_pheno(c,value) result(res) 1132 integer ,intent(in) :: c 1133 character(len=LEN_DEF) ,intent(in) :: value 1134 integer(kind=KIND_PHENO) :: res 1135 integer :: nbvalue 1136 integer :: i,ch 1137 character(len=LEN_DEF) , dimension (:,:) ,allocatable :: buf_value_pheno 1138 res = VAL_MIN_INDEX_PHENO 1139 1140 do i=1,nb_value_pheno(c) 1141 if ( value_pheno(c,i) == trim(value) ) then 1142 return 1143 end if 1144 res = res + 1 1145 end do 1146 1147 !manage error 1148 if ( i > VAL_MAX_INDEX_PHENO ) then 1149 nbvalue = VAL_MAX_INDEX_PHENO 1150 nbvalue = nbvalue - VAL_MIN_INDEX_PHENO 1151 1152 write (0,*) "* Too many allele detected for internal data structure " 1153 write (0,*) "number max [",nbvalue,']' 1154 do ch=1,nchr 1155 print *,"------chromo : "//trim(chromo(ch))//"---------- size: ",nb_value_pheno(ch) 1156 print *,(trim(value_pheno(ch,i))//",",i=1,nb_value_pheno(ch)) 1157 end do 1158 stop 1159 end if 1160 1161 if ( size(value_pheno,2) <= nb_value_pheno(c) ) then 1162 allocate (buf_value_pheno(nchr,size(value_pheno,2))) 1163 buf_value_pheno = value_pheno 1164 deallocate(value_pheno) 1165 allocate( value_pheno(nchr,size(buf_value_pheno,2)+256) ) 1166 value_pheno="" 1167 do ch=1,nchr 1168 do i=1,nb_value_pheno(ch) 1169 value_pheno(ch,i)=trim(buf_value_pheno(ch,i)) 1170 end do 1171 end do 1172 deallocate(buf_value_pheno) 1173 end if 1174 1175 nb_value_pheno(c)=nb_value_pheno(c)+1 1176 value_pheno(c,nb_value_pheno(c))=trim(value) 1177 return 1178 1179 end function set_pheno
stop_application
[ Top ] [ m_qtlmap_tools ] [ Subroutines ]
NAME
stop_application
DESCRIPTION
stop the application with an error message
INPUTS
message : an error message
SOURCE
224 subroutine stop_application(message) 225 character(len=*),intent(in) ::message 226 call log_mess('*** message error :'//message//' ***',ERROR_DEF) 227 stop 'runtime error' 228 229 end subroutine stop_application
stop_on_error
[ Top ] [ m_qtlmap_tools ] [ Subroutines ]
NAME
stop_on_error
DESCRIPTION
stop the application with an error message on a read file context
INPUTS
iostat : iostat error code file : the pathname file lineNumber : the number of the current line message : a message error
SOURCE
199 subroutine stop_on_error(iostat,file,lineNumber,message) 200 integer , intent(in) ::iostat 201 character(len=*),intent(in) ::file 202 integer , intent(in) ::lineNumber 203 character(len=*),intent(in) ::message 204 205 if (iostat /= 0 .and. iostat/=-2) then 206 call log_mess('*** Read error ['//trim(file)//'] ' // & 207 'line '//trim(iget_str(lineNumber))//':bad specification [message:'//adjustr(message)//'] ***',ERROR_DEF) 208 stop 'read error' 209 endif 210 211 end subroutine stop_on_error
str
[ Top ] [ m_qtlmap_tools ] [ Subroutines ]
NAME
str
DESCRIPTION
generic interface to convert a float/double or an integer to a string
SOURCE
123 interface str 124 module procedure iget_str, & 125 fget_str, & 126 dget_str 127 end interface str
xaldane
[ Top ] [ m_qtlmap_tools ] [ Subroutines ]
NAME
xaldane
DESCRIPTION
INPUTS
d :
RETURN
NOTES
Calcul du taux de recombinaison en fonction de la distance de xaldane
SOURCE
599 function xaldane(do) result (val) 600 real (kind=dp) , intent(in) :: do 601 real (kind=dp) :: val 602 val=0.5d0*(1.d0-dexp(-2.d0*do)) 603 return 604 end function xaldane
xosambi
[ Top ] [ m_qtlmap_tools ] [ Subroutines ]
NAME
xosambi
DESCRIPTION
INPUTS
d :
RETURN
NOTES
Calcul du taux de recombinaison en fonction de la distance de Kosambi
SOURCE
621 function xosambi(do) result (val) 622 real (kind=dp) , intent(in) :: do 623 real (kind=dp) :: val,z 624 625 z=dexp(4.d0*do) 626 val=0.5d0*((z-1.d0)/(z+1.d0)) 627 return 628 end function xosambi