JHUGen MELA  v2.4.1
Matrix element calculations as used in JHUGen. MELA is an important tool that was used for the Higgs boson discovery and for precise measurements of its structure and interactions. Please see the website https://spin.pha.jhu.edu/ and papers cited there for more details, and kindly cite those papers when using this code.
coli_aux2.F90
Go to the documentation of this file.
1 !!
2 !! File coli_aux2.F90 is part of COLLIER
3 !! - A Complex One-Loop Library In Extended Regularizations
4 !!
5 !! Copyright (C) 2015, 2016 Ansgar Denner, Stefan Dittmaier, Lars Hofer
6 !!
7 !! COLLIER is licenced under the GNU GPL version 3, see COPYING for details.
8 !!
9 
10 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
11 !
12 ! ***********************
13 ! * module coli_aux2 *
14 ! * by Lars Hofer *
15 ! ***********************
16 !
17 ! functions and subroutines:
18 !
19 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
20 
21 #include "global_coli.h"
22 
23 module coli_aux2
24 
25  use master
26 
27  implicit none
28 
31  integer :: maxerrout_coli=100
32  integer :: stdout_coli=6,mode_coli=0,closed_coli=-999
33  double precision :: dprec_coli, reqacc_coli, critacc_coli
34 
35  double precision :: acc_inf, acc_def_b, acc_def_c0, acc_def_d0
36  double precision :: acc_req_c, acc_req_d, acc_req_cind
37  double precision :: impest_d=1d1, impest_c=1d1
38  integer :: rmax_b, rmax_c, rmax_d
39  logical :: inflev_coli
40 
41 contains
42 
43 
44  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
45  ! subroutine initcoli_in_collier
46  ! fixing of default values for various variables *
47  !
48  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
49 
50  subroutine initcoli_in_collier
51 
52  implicit none
53 ! logical init
54 !
55 ! data init /.false./
56 ! save init
57 !
58 ! if (init) return
59 ! init = .true.
60 
61  ! scale factor for *small* masses masses
62  call setminfscale2_coli(1d0)
63 #ifdef SING
64  ! shift of mass-singular squared logarithms
65  call setshiftms2_coli(0d0)
66 #endif
67 
68 
69  ! estimate of CPU precision
70  ! calacc = dprec_coli
71  ! infinitesimal parameter
72  ! eps = dprec_coli/4d0
73  ! size of imaginary parts below which explicit ieps take over
74  ! impacc = dprec_coli/4d4
75  ! done via GetCPUprec_cll
76  ! call setprecpars_coli(dprec_coli)
77 
78  ! infinitesimal accuracy
79  acc_inf = 1d40
80 
81 ! acc_min_C = 1d0
82 ! acc_min_D = 1d0
83 
84 
85  end subroutine initcoli_in_collier
86 
87 
88 
89  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
90  ! subroutine InitErrCnt_coli(val)
91  !
92  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
93 
94  subroutine initerrcnt_coli(val)
95 
96  integer, intent(in) :: val
97 
98  errcnt_coli = val
99 
100  end subroutine initerrcnt_coli
101 
102 
103 
104 
105 
106  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
107  ! subroutine InitCritPointsCnt_coli(val)
108  !
109  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
110 
111  subroutine initcritpointscnt_coli(val)
112 
113  integer, intent(in) :: val
114 
115  critpointscnt_coli = val
116 
117  end subroutine initcritpointscnt_coli
118 
119 
120 
121 
122 
123  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
124  ! subroutine SetErrFlag_coli(err)
125  !
126  ! -1 Check failed
127  ! -4 Argument on cut (only if CHECK set)
128  ! -5 Critical event
129  ! -6 No reduction method works
130  ! -6 momenta not 4-dimensional
131  ! -7 specific numerical problem
132  !
133  ! -10 Case not supported/implemented
134  !
135  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
136 
137  subroutine seterrflag_coli(err)
138 
139  integer, intent(in) :: err
140 
141  errflag_coli = err
142 
143  end subroutine seterrflag_coli
144 
145 
146 
147 
148 
149  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
150  ! subroutine getErrFlag_coli(err)
151  !
152  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
153 
154  subroutine geterrflag_coli(err)
155 
156  integer, intent(out) :: err
157 
158  err = errflag_coli
159 
160  end subroutine geterrflag_coli
161 
162 
163 
164 
165 
166  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
167  ! subroutine setnerrout_coli(nerrout)
168  !
169  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
170 
171  subroutine setnerrout_coli(nerrout)
172 
173  integer, intent(in) :: nerrout
174 
175  nerrout_coli = nerrout
176 
177  if(inflev_coli.and.ninfout_coli.ne.closed_coli)then
178  write(ninfout_coli,*) 'COLI: nerrout_coli set to = ',nerrout_coli
179  endif
180 
181  end subroutine setnerrout_coli
182 
183 
184 
185 
186 
187  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
188  ! subroutine setncpout_coli(ncpout)
189  !
190  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
191 
192  subroutine setncpout_coli(ncpout)
193 
194  integer, intent(in) :: ncpout
195 
196  ncpout_coli = ncpout
197 
198  if(inflev_coli.and.ninfout_coli.ne.closed_coli)then
199  write(ninfout_coli,*) 'COLI: ncpout_coli set to = ',ncpout_coli
200  endif
201 
202  end subroutine setncpout_coli
203 
204 
205 
206 
207 
208  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
209  ! subroutine setninfout_coli(ninfout)
210  !
211  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
212 
213  subroutine setninfout_coli(ninfout)
214 
215  integer, intent(in) :: ninfout
216 
217  ninfout_coli = ninfout
218 
219  if(inflev_coli.and.ninfout_coli.ne.closed_coli)then
220  write(ninfout_coli,*) 'COLI: ninfout_coli set to = ',ninfout_coli
221  endif
222 
223  end subroutine setninfout_coli
224 
225 
226 
227 
228 
229  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
230  ! subroutine setnstatsout_coli(nstatsout)
231  !
232  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
233 
234  subroutine setnstatsout_coli(nstatsout)
235 
236  integer, intent(in) :: nstatsout
237 
238  nstatsout_coli = nstatsout
239 
240  if(inflev_coli.and.ninfout_coli.ne.closed_coli)then
241  write(ninfout_coli,*) 'COLI: nstatsout_coli set to = ',nstatsout_coli
242  endif
243 
244  end subroutine setnstatsout_coli
245 
246 
247 
248 
249 
250  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
251  ! subroutine SetMaxErrOut_coli(errmax)
252  !
253  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
254 
255  subroutine setmaxerrout_coli(errmax)
256 
257  integer, intent(in) :: errmax
258 
259  maxerrout_coli = errmax
260 
261  if(inflev_coli.and.ninfout_coli.ne.closed_coli)then
262  write(ninfout_coli,*) 'COLI: maxErrOut_coli set to = ',maxerrout_coli
263  endif
264 
265  end subroutine setmaxerrout_coli
266 
267 
268 
269 
270 
271  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
272  ! subroutine SetErroutlev_coli(erroutlev)
273  !
274  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
275 
276  subroutine seterroutlev_coli(erroutlev)
277 
278  integer, intent(in) :: erroutlev
279 
280  erroutlev_coli = erroutlev
281 
282 
283  end subroutine seterroutlev_coli
284 
285 
286 
287 
288 
289  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
290  ! subroutine ErrOut_coli(sub,err,flag)
291  !
292  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
293 
294 ! Suppression of output must be implemented in calling routines!
295 
296  subroutine errout_coli(sub,err,flag)
297 
298  character(len=*), intent(in) :: sub, err
299 ! integer, parameter :: maxErrOut=100
300  logical, intent(out) :: flag
301 
302  flag = .false.
303  if(erroutlev_coli.eq.0) return
304 
305  errcnt_coli = errcnt_coli + 1
306  if (errcnt_coli.le.maxerrout_coli) then
307  write(nerrout_coli,*)
308  write(nerrout_coli,*)
309  write(nerrout_coli,*)
310  write(nerrout_coli,*) '***********************************************************'
311  write(nerrout_coli,*) 'ERROR NO.', errcnt_coli
312  write(nerrout_coli,*) 'in routine: ', trim(sub)
313  write(nerrout_coli,*) trim(err)
314  call writemaster_cll(nerrout_coli)
315  flag = .true.
316  elseif (errcnt_coli.eq.maxerrout_coli+1) then
317  write(nerrout_coli,*)
318  write(nerrout_coli,*)
319  write(nerrout_coli,*)
320  write(nerrout_coli,*) '***********************************************************'
321  write(nerrout_coli,*)
322  write(nerrout_coli,*) ' Further output of Errors will be suppressed '
323  write(nerrout_coli,*)
324  endif
325 
326  end subroutine errout_coli
327 
328 
329 
330 
331 
332  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
333  ! subroutine CritPointsOut_coli(sub,acc)
334  !
335  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
336 
337  subroutine critpointsout_coli(sub,acc)
338 
339  character(len=*), intent(in) :: sub
340  double precision, intent(in) :: acc
341  integer :: i
342 
343 #include "common_coli.h"
344 
346 
347  write(ncpout_coli,*)
348  write(ncpout_coli,*)
349  write(ncpout_coli,*)
350  write(ncpout_coli,*) '***********************************************************'
351  write(ncpout_coli,*) 'Critical Point NO.', critpointscnt_coli
352  write(ncpout_coli,*) 'in integral: ', trim(sub)
353  write(ncpout_coli,*) 'estimated accuracy: ', acc
354  write(ncpout_coli,*) '-----------------------------------------------------------'
355  write(ncpout_coli,*) 'GLOBAL PARAMETERS:'
356  write(ncpout_coli,*) 'muUV2 = ', muuv2
357  write(ncpout_coli,*) 'muIR2 = ', muir2
358 #ifdef SING
359  write(ncpout_coli,*) 'deltaUV = ', deltauv
360  write(ncpout_coli,*) 'deltaIR1 = ', delta1ir
361  write(ncpout_coli,*) 'deltaIR2 = ', delta2ir
362 #endif
363  write(ncpout_coli,*) 'nminf = ', ncoliminf
364  do i=1,ncoliminf
365  write(ncpout_coli,*) 'minf2 = ', i, coliminf2(i)
366  end do
367  write(ncpout_coli,*) 'dprec = ', dprec_coli
368  write(ncpout_coli,*) 'reqacc = ', reqacc_coli
369  write(ncpout_coli,*) 'critacc = ', critacc_coli
370  write(ncpout_coli,*) 'ErrFlag = ', errflag_coli
371 ! write(ncpout_coli,*) '------------------------------------------------------------'
372  call writemaster_cll(ncpout_coli)
373 
374  end subroutine critpointsout_coli
375 
376 
377 
378 
379 
380 
381 
382 
383  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
384  ! subroutine SetMode_coli(mode)
385  !
386  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
387 
388  subroutine setmode_coli(mode)
389 
390  integer, intent(in) :: mode
391 
392  mode_coli = mode
393  if(mode.lt.1) then
394  impest_c = 1d1
395  impest_d = 1d1
396  else
397  impest_c = 1d0
398  impest_d = 1d0
399  end if
400 
401  if(inflev_coli.and.ninfout_coli.ne.closed_coli)then
402  write(ninfout_coli,*) 'COLI: Mode_coli set to = ',mode_coli
403  endif
404 
405  end subroutine setmode_coli
406 
407 
408 
409 
410 
411  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
412  ! subroutine setprec_coli(dprec)
413  !
414  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
415 
416  subroutine setprec_coli(dprec)
417 
418  double precision, intent(in) :: dprec
419 
420  dprec_coli = dprec
421 
422 ! adapted to precision of D0 and C0
424  acc_def_c0 = 1d1*dprec_coli
425  acc_def_d0 = 1d1*dprec_coli
426 
427  call setprecpars_coli(dprec)
428 
429  if(inflev_coli.and.ninfout_coli.ne.closed_coli)then
430  write(ninfout_coli,*) 'COLI: dprec_coli set to = ',dprec_coli
431  endif
432 
433  end subroutine setprec_coli
434 
435 
436 
437 
438 
439  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
440  ! subroutine SetAcc_coli(reqacc,critacc)
441  !
442  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
443 
444  subroutine setacc_coli(reqacc,critacc)
445 
446  double precision, intent(in) :: reqacc,critacc
447 
448  reqacc_coli = reqacc
449  acc_req_c = reqacc
450  acc_req_cind = reqacc/1d1
451  acc_req_d = reqacc
452 
453  critacc_coli = critacc
454 
455  if(inflev_coli.and.ninfout_coli.ne.closed_coli)then
456  write(ninfout_coli,*) 'COLI: reqacc_coli set to = ',reqacc_coli
457  write(ninfout_coli,*) 'COLI: critacc_coli set to = ',critacc_coli
458  endif
459 
460  end subroutine setacc_coli
461 
462 
463 
464 
465 
466  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
467  ! subroutine SetReqAcc_coli(reqacc)
468  !
469  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
470 
471  subroutine setreqacc_coli(reqacc)
472 
473  double precision, intent(in) :: reqacc
474 
475  reqacc_coli = reqacc
476  acc_req_c = reqacc
477  acc_req_cind = reqacc/1d1
478  acc_req_d = reqacc
479 
480  if(inflev_coli.and.ninfout_coli.ne.closed_coli)then
481  write(ninfout_coli,*) 'COLI: reqacc_coli set to = ',reqacc_coli
482  endif
483 
484  end subroutine setreqacc_coli
485 
486 
487 
488 
489 
490  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
491  ! subroutine SetCritAcc_coli(critacc)
492  !
493  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
494 
495  subroutine setcritacc_coli(critacc)
496 
497  double precision, intent(in) :: critacc
498 
499  critacc_coli = critacc
500 
501  if(inflev_coli.and.ninfout_coli.ne.closed_coli)then
502  write(ninfout_coli,*) 'COLI: critacc_coli set to = ',critacc_coli
503  endif
504 
505  end subroutine setcritacc_coli
506 
507 
508 
509 
510 
511  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
512  ! subroutine SetRitmax_coli(ritmax_B,ritmax_C,ritmax_D)
513  !
514  !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
515 
516  subroutine setritmax_coli(ritmax_B,ritmax_C,ritmax_D)
517 
518  integer, intent(in) :: ritmax_B,ritmax_C, ritmax_D
519 
520  if (ritmax_d.lt.6) then
521  if(inflev_coli.and.ninfout_coli.ne.closed_coli)then
522  write(ninfout_coli,*) 'COLI: SetRitmax_coli'
523  write(ninfout_coli,*) 'rmax_D has to be at least 6 --> it is set to 6'
524  end if
525  rmax_d = 6
526  else
527  rmax_d = ritmax_d
528  if(inflev_coli.and.ninfout_coli.ne.closed_coli)then
529  write(ninfout_coli,*) 'COLI: rmax_D set to = ',rmax_d
530  endif
531  end if
532 
533  if (ritmax_c.le.rmax_d) then
534  if(inflev_coli.and.ninfout_coli.ne.closed_coli)then
535  write(ninfout_coli,*) 'COLI: SetRitmax_coli'
536  write(ninfout_coli,*) 'rmax_C has to be larger than rmax_D --> it is set to rmax_D+1'
537  end if
538  rmax_c = rmax_d+1
539  else
540  rmax_c = ritmax_c
541  if(inflev_coli.and.ninfout_coli.ne.closed_cll)then
542  write(ninfout_coli,*) 'COLI: rmax_C set to = ',rmax_c
543  endif
544  end if
545 
546  if (ritmax_b.le.rmax_c) then
547  if(inflev_coli.and.ninfout_coli.ne.closed_coli)then
548  write(ninfout_coli,*) 'COLI: SetRitmax_coli'
549  write(ninfout_coli,*) 'rmax_B has to be larger than rmax_C --> it is set to rmax_C+1'
550  end if
551  rmax_b = rmax_c+1
552  else
553  rmax_b = ritmax_b
554  if(inflev_coli.and.ninfout_coli.ne.closed_coli)then
555  write(ninfout_coli,*) 'COLI: rmax_B set to = ',rmax_b
556  endif
557  end if
558 
559 
560  end subroutine setritmax_coli
561 
562 
563 end module coli_aux2
coli_aux2::initcritpointscnt_coli
subroutine initcritpointscnt_coli(val)
Definition: coli_aux2.F90:112
coli_aux2::acc_def_d0
double precision acc_def_d0
Definition: coli_aux2.F90:35
coli_aux2::setncpout_coli
subroutine setncpout_coli(ncpout)
Definition: coli_aux2.F90:193
coli_aux2::setcritacc_coli
subroutine setcritacc_coli(critacc)
Definition: coli_aux2.F90:496
coli_aux2::nerrout_coli
integer nerrout_coli
Definition: coli_aux2.F90:29
coli_aux2::rmax_d
integer rmax_d
Definition: coli_aux2.F90:38
coli_aux2::initerrcnt_coli
subroutine initerrcnt_coli(val)
Definition: coli_aux2.F90:95
coli_aux2::acc_req_c
double precision acc_req_c
Definition: coli_aux2.F90:36
coli_aux2::dprec_coli
double precision dprec_coli
Definition: coli_aux2.F90:33
coli_aux2::setnstatsout_coli
subroutine setnstatsout_coli(nstatsout)
Definition: coli_aux2.F90:235
coli_aux2::stdout_coli
integer stdout_coli
Definition: coli_aux2.F90:32
coli_aux2::reqacc_coli
double precision reqacc_coli
Definition: coli_aux2.F90:33
coli_aux2::initcoli_in_collier
subroutine initcoli_in_collier
Definition: coli_aux2.F90:51
coli_aux2::setmaxerrout_coli
subroutine setmaxerrout_coli(errmax)
Definition: coli_aux2.F90:256
coli_aux2
Definition: coli_aux2.F90:23
coli_aux2::ninfout_coli
integer ninfout_coli
Definition: coli_aux2.F90:29
coli_aux2::nstatsout_coli
integer nstatsout_coli
Definition: coli_aux2.F90:29
coli_aux2::acc_inf
double precision acc_inf
Definition: coli_aux2.F90:35
coli_aux2::setritmax_coli
subroutine setritmax_coli(ritmax_B, ritmax_C, ritmax_D)
Definition: coli_aux2.F90:517
coli_aux2::critacc_coli
double precision critacc_coli
Definition: coli_aux2.F90:33
coli_aux2::geterrflag_coli
subroutine geterrflag_coli(err)
Definition: coli_aux2.F90:155
coli_aux2::impest_c
double precision impest_c
Definition: coli_aux2.F90:37
coli_aux2::acc_def_c0
double precision acc_def_c0
Definition: coli_aux2.F90:35
coli_aux2::errflag_coli
integer errflag_coli
Definition: coli_aux2.F90:30
coli_aux2::impest_d
double precision impest_d
Definition: coli_aux2.F90:37
coli_aux2::setreqacc_coli
subroutine setreqacc_coli(reqacc)
Definition: coli_aux2.F90:472
coli_aux2::errcnt_coli
integer errcnt_coli
Definition: coli_aux2.F90:30
coli_aux2::acc_req_d
double precision acc_req_d
Definition: coli_aux2.F90:36
coli_aux2::setnerrout_coli
subroutine setnerrout_coli(nerrout)
Definition: coli_aux2.F90:172
coli_aux2::acc_req_cind
double precision acc_req_cind
Definition: coli_aux2.F90:36
coli_aux2::setninfout_coli
subroutine setninfout_coli(ninfout)
Definition: coli_aux2.F90:214
coli_aux2::mode_coli
integer mode_coli
Definition: coli_aux2.F90:32
coli_aux2::critpointscnt_coli
integer critpointscnt_coli
Definition: coli_aux2.F90:30
coli_aux2::erroutlev_coli
integer erroutlev_coli
Definition: coli_aux2.F90:30
coli_aux2::ncpout_coli
integer ncpout_coli
Definition: coli_aux2.F90:29
coli_aux2::maxerrout_coli
integer maxerrout_coli
Definition: coli_aux2.F90:31
coli_aux2::critpointsout_coli
subroutine critpointsout_coli(sub, acc)
Definition: coli_aux2.F90:338
coli_aux2::setprec_coli
subroutine setprec_coli(dprec)
Definition: coli_aux2.F90:417
coli_aux2::setmode_coli
subroutine setmode_coli(mode)
Definition: coli_aux2.F90:389
coli_aux2::seterroutlev_coli
subroutine seterroutlev_coli(erroutlev)
Definition: coli_aux2.F90:277
coli_aux2::closed_coli
integer closed_coli
Definition: coli_aux2.F90:32
coli_aux2::rmax_c
integer rmax_c
Definition: coli_aux2.F90:38
coli_aux2::inflev_coli
logical inflev_coli
Definition: coli_aux2.F90:39
coli_aux2::setacc_coli
subroutine setacc_coli(reqacc, critacc)
Definition: coli_aux2.F90:445
coli_aux2::errout_coli
subroutine errout_coli(sub, err, flag)
Definition: coli_aux2.F90:297
coli_aux2::seterrflag_coli
subroutine seterrflag_coli(err)
Definition: coli_aux2.F90:138
coli_aux2::acc_def_b
double precision acc_def_b
Definition: coli_aux2.F90:35
coli_aux2::rmax_b
integer rmax_b
Definition: coli_aux2.F90:38