lookup_table.F90 coverage: 83.33 %func 72.44 %block
1) module Lookup_Table_module
2)
3) use PFLOTRAN_Constants_module
4)
5) implicit none
6)
7) private
8)
9) #include "petsc/finclude/petscsys.h"
10)
11) type, abstract, public :: lookup_table_base_type
12) PetscInt :: dim
13) PetscInt :: dims(3)
14) PetscReal, pointer :: data(:)
15) class(lookup_table_axis_type), pointer :: axis1
16) contains
17) procedure(LookupTableEvaluateDummy), deferred, public :: Sample
18) end type lookup_table_base_type
19)
20) type, public, extends(lookup_table_base_type) :: lookup_table_uniform_type
21) class(lookup_table_axis_type), pointer :: axis2
22) class(lookup_table_axis_type), pointer :: axis3
23) contains
24) procedure, public :: Sample => LookupTableEvaluateUniform
25) end type lookup_table_uniform_type
26)
27) type, public, extends(lookup_table_base_type) :: lookup_table_general_type
28) class(lookup_table_axis2_general_type), pointer :: axis2
29) contains
30) procedure, public :: Sample => LookupTableEvaluateGeneral
31) end type lookup_table_general_type
32)
33) type, public :: lookup_table_axis_type
34) PetscInt :: itype
35) PetscInt :: saved_index
36) PetscReal, pointer :: values(:)
37) end type lookup_table_axis_type
38)
39) type, public, extends(lookup_table_axis_type) :: lookup_table_axis2_general_type
40) PetscInt :: saved_index2
41) end type lookup_table_axis2_general_type
42)
43) abstract interface
44) function LookupTableEvaluateDummy(this,lookup1,lookup2,lookup3)
45) import lookup_table_base_type
46) implicit none
47) class(lookup_table_base_type) :: this
48) PetscReal :: lookup1
49) PetscReal, optional :: lookup2
50) PetscReal, optional :: lookup3
51) PetscReal :: LookupTableEvaluateDummy
52) end function LookupTableEvaluateDummy
53) end interface
54)
55) interface LookupTableTest
56) module procedure LookupTableTest1D
57) module procedure LookupTableTest2D
58) end interface
59)
60) interface LookupTableDestroy
61) module procedure LookupTableUniformDestroy
62) module procedure LookupTableGeneralDestroy
63) end interface
64)
65) public :: LookupTableCreateUniform, &
66) LookupTableCreateGeneral, &
67) LookupTableDestroy, &
68) LookupTableTest
69)
70) contains
71)
72) ! ************************************************************************** !
73)
74) subroutine LookupTableBaseInit(lookup_table)
75) !
76) ! Author: Glenn Hammond
77) ! Date: 10/15/14
78) !
79)
80) implicit none
81)
82) class(lookup_table_base_type) :: lookup_table
83)
84) lookup_table%dim = 0
85) lookup_table%dims = 0
86) nullify(lookup_table%data)
87) nullify(lookup_table%axis1)
88)
89) end subroutine LookupTableBaseInit
90)
91) ! ************************************************************************** !
92)
93) function LookupTableCreateUniform(dim)
94) !
95) ! Author: Glenn Hammond
96) ! Date: 10/15/14
97) !
98)
99) implicit none
100)
101) PetscInt :: dim
102)
103) class(lookup_table_uniform_type), pointer :: LookupTableCreateUniform
104)
105) class(lookup_table_uniform_type), pointer :: lookup_table
106)
107) allocate(lookup_table)
108) call LookupTableBaseInit(lookup_table)
109) lookup_table%dim = dim
110) nullify(lookup_table%axis2)
111) nullify(lookup_table%axis3)
112) allocate(lookup_table%axis1)
113) call LookupTableAxisInit(lookup_table%axis1)
114) if (dim > 1) then
115) allocate(lookup_table%axis2)
116) call LookupTableAxisInit(lookup_table%axis2)
117) endif
118) if (dim > 2) then
119) allocate(lookup_table%axis3)
120) call LookupTableAxisInit(lookup_table%axis3)
121) endif
122)
123) LookupTableCreateUniform => lookup_table
124)
125) end function LookupTableCreateUniform
126)
127) ! ************************************************************************** !
128)
129) function LookupTableCreateGeneral(dim)
130) !
131) ! Author: Glenn Hammond
132) ! Date: 10/15/14
133) !
134)
135) implicit none
136)
137) PetscInt :: dim
138)
139) class(lookup_table_general_type), pointer :: LookupTableCreateGeneral
140)
141) class(lookup_table_general_type), pointer :: lookup_table
142)
143) allocate(lookup_table)
144) call LookupTableBaseInit(lookup_table)
145) lookup_table%dim = dim
146) nullify(lookup_table%axis2)
147) allocate(lookup_table%axis1)
148) call LookupTableAxisInit(lookup_table%axis1)
149) if (dim > 1) then
150) allocate(lookup_table%axis2)
151) call LookupTableAxisInit(lookup_table%axis2)
152) lookup_table%axis2%saved_index2 = 1
153) endif
154)
155) LookupTableCreateGeneral => lookup_table
156)
157) end function LookupTableCreateGeneral
158)
159) ! ************************************************************************** !
160)
161) subroutine LookupTableAxisInit(axis)
162) !
163) ! Author: Glenn Hammond
164) ! Date: 10/15/14
165) !
166)
167) implicit none
168)
169) class(lookup_table_axis_type) :: axis
170)
171) axis%itype = 0
172) axis%saved_index = 1
173) nullify(axis%values)
174)
175) end subroutine LookupTableAxisInit
176)
177) ! ************************************************************************** !
178)
179) function LookupTableEvaluateUniform(this,lookup1,lookup2,lookup3)
180) !
181) ! Author: Glenn Hammond
182) ! Date: 10/15/14
183) !
184) implicit none
185)
186) class(lookup_table_uniform_type) :: this
187) PetscReal :: lookup1
188) PetscReal, optional :: lookup2
189) PetscReal, optional :: lookup3
190)
191) PetscReal :: LookupTableEvaluateUniform
192)
193) call LookupTableIndexUniform(this,lookup1,lookup2,lookup3)
194) if (present(lookup3)) then
195) ! call LookupTableInterpolate3DUniform(this,lookup1,lookup2,lookup3,LookupTableEvaluateUniform)
196) else if (present(lookup2)) then
197) call LookupTableInterpolate2DUniform(this,lookup1,lookup2,LookupTableEvaluateUniform)
198) else
199) call LookupTableInterpolate1D(this,lookup1,LookupTableEvaluateUniform)
200) endif
201)
202) end function LookupTableEvaluateUniform
203)
204)
205) ! ************************************************************************** !
206)
207) function LookupTableEvaluateGeneral(this,lookup1,lookup2,lookup3)
208) !
209) ! Author: Glenn Hammond
210) ! Date: 10/15/14
211) !
212) implicit none
213)
214) class(lookup_table_general_type) :: this
215) PetscReal :: lookup1
216) PetscReal, optional :: lookup2
217) PetscReal, optional :: lookup3
218)
219) PetscReal :: LookupTableEvaluateGeneral
220)
221) call LookupTableIndexGeneral(this,lookup1,lookup2,lookup3)
222) if (present(lookup3)) then
223) ! call LookupTableInterpolate3DGeneral(this,lookup1,lookup2,lookup3,LookupTableEvaluateGeneral)
224) else if (present(lookup2)) then
225) call LookupTableInterpolate2DGeneral(this,lookup1,lookup2,LookupTableEvaluateGeneral)
226) else
227) call LookupTableInterpolate1D(this,lookup1,LookupTableEvaluateGeneral)
228) endif
229)
230) end function LookupTableEvaluateGeneral
231)
232) ! ************************************************************************** !
233)
234) subroutine LookupTableIndexUniform(this,lookup1,lookup2,lookup3)
235) !
236) ! Author: Glenn Hammond
237) ! Date: 10/15/14
238) !
239) implicit none
240)
241) class(lookup_table_uniform_type) :: this
242) PetscReal :: lookup1
243) PetscReal, optional :: lookup2
244) PetscReal, optional :: lookup3
245)
246) call LookupTableAxisIndexUniform(this%axis1,lookup1)
247) if (associated(this%axis2)) then
248) call LookupTableAxisIndexUniform(this%axis2,lookup2)
249) endif
250) if (associated(this%axis3)) then
251) call LookupTableAxisIndexUniform(this%axis3,lookup3)
252) endif
253)
254) end subroutine LookupTableIndexUniform
255)
256) ! ************************************************************************** !
257)
258) subroutine LookupTableIndexGeneral(this,lookup1,lookup2,lookup3)
259) !
260) ! Author: Glenn Hammond
261) ! Date: 10/15/14
262) !
263) implicit none
264)
265) class(lookup_table_general_type) :: this
266) PetscReal :: lookup1
267) PetscReal, optional :: lookup2
268) PetscReal, optional :: lookup3
269)
270) PetscInt :: ja, jb
271) PetscInt :: istart, iend
272) class(lookup_table_axis2_general_type), pointer :: axis2
273)
274) ! axis 1 corresponds to the j dim below
275) call LookupTableAxisIndexGeneral(lookup1,this%axis1%values, &
276) this%axis1%saved_index)
277)
278) if (.not. associated(this%axis2)) return
279)
280) axis2 => this%axis2
281)
282) ja = this%axis1%saved_index
283) if (ja > 0) then
284) jb = max(min(ja+1,this%dims(1)),1)
285) else
286) ja = 1
287) jb = 1
288) endif
289)
290) iend = ja*this%dims(2)
291) istart = iend - this%dims(2) + 1
292) call LookupTableAxisIndexGeneral(lookup2,axis2%values(istart:iend), &
293) axis2%saved_index)
294) if (ja /= jb) then
295) iend = jb*this%dims(2)
296) istart = iend - this%dims(2) + 1
297) call LookupTableAxisIndexGeneral(lookup2,axis2%values(istart:iend), &
298) axis2%saved_index2)
299) else
300) axis2%saved_index2 = axis2%saved_index
301) endif
302)
303) end subroutine LookupTableIndexGeneral
304)
305) ! ************************************************************************** !
306)
307) subroutine LookupTableAxisIndexUniform(this,lookup1)
308) !
309) ! Author: Glenn Hammond
310) ! Date: 10/15/14
311) !
312) implicit none
313)
314) class(lookup_table_axis_type) :: this
315) PetscReal :: lookup1
316)
317) PetscInt :: i1
318) PetscInt :: size1
319) PetscReal :: begin1
320)
321) size1 = size(this%values)
322) begin1 = this%values(1)
323) i1 = int((lookup1 - begin1) / (this%values(size1) - begin1) * (size1-1) + 1)
324) ! truncate i1 to zero indicating the value is below the range specified
325) i1 = max(min(i1,size1),0)
326) this%saved_index = i1
327)
328) end subroutine LookupTableAxisIndexUniform
329)
330) ! ************************************************************************** !
331)
332) subroutine LookupTableAxisIndexGeneral(lookup1,values,saved_index)
333) !
334) ! Author: Glenn Hammond
335) ! Date: 10/15/14
336) !
337) implicit none
338)
339) PetscReal :: lookup1
340) PetscInt :: saved_index
341) PetscReal :: values(:)
342)
343) PetscInt :: i1
344) PetscInt :: j1
345) PetscInt :: mid1
346) PetscInt :: sizei
347)
348) sizei = size(values)
349) if (lookup1 < values(1)) then
350) saved_index = 0
351) return
352) else if (lookup1 > values(sizei)) then
353) saved_index = sizei
354) return
355) endif
356) i1 = max(min(saved_index,sizei-1),1)
357) if (lookup1 > values(i1+1) .or. &
358) lookup1 < values(i1)) then
359) ! move either up or down array
360) if (lookup1 > values(i1+1)) then
361) i1 = i1+1
362) j1 = sizei
363) else
364) j1 = i1
365) i1 = 1
366) endif
367) do
368) mid1 = (j1+i1) / 2
369) if (lookup1 > values(mid1)) then
370) i1 = mid1
371) else
372) j1 = mid1
373) endif
374) if (j1-i1 <= 1) exit
375) enddo
376) endif
377) saved_index = i1
378)
379) end subroutine LookupTableAxisIndexGeneral
380)
381) ! ************************************************************************** !
382)
383) subroutine LookupTableInterpolate1D(this,lookup1,result)
384) !
385) ! Author: Glenn Hammond
386) ! Date: 10/15/14
387) !
388) use Utility_module, only : Interpolate
389)
390) implicit none
391)
392) class(lookup_table_base_type) :: this
393) PetscReal :: lookup1
394) PetscReal :: result
395)
396) PetscInt :: i1, j1
397)
398) i1 = this%axis1%saved_index
399) if (i1 > 0) then
400) j1 = max(min(i1+1,this%dims(1)),1)
401) call Interpolate(this%axis1%values(j1),this%axis1%values(i1),lookup1, &
402) this%data(j1),this%data(i1),result)
403) else ! catch values below axis range
404) result = this%data(1)
405) endif
406)
407) end subroutine LookupTableInterpolate1D
408)
409) ! ************************************************************************** !
410)
411) subroutine LookupTableInterpolate2DUniform(this,lookup1,lookup2,result)
412) !
413) ! Author: Glenn Hammond
414) ! Date: 10/15/14
415) !
416) use Utility_module, only : Interpolate, InterpolateBilinear
417)
418) implicit none
419)
420) class(lookup_table_uniform_type) :: this
421) PetscReal :: lookup1
422) PetscReal :: lookup2
423) PetscReal :: result
424)
425) PetscInt :: i1, i2, j1, j2
426) PetscReal :: x1, x2, y1, y2, z1, z2, z3, z4
427) PetscInt :: sizei, sizej
428)
429) ! x1,y2,z3 ------ x2,y2,z4
430) ! | |
431) ! | |
432) ! | x,y |
433) ! | |
434) ! x1,y1,z1 ------ x2,y1,z2
435)
436) result = UNINITIALIZED_DOUBLE
437) sizei = this%dims(1)
438) sizej = this%dims(2)
439) i1 = this%axis1%saved_index
440) j1 = this%axis2%saved_index
441) ! index axes
442) if (i1 > 0) then
443) i2 = max(min(i1+1,sizei),1)
444) else
445) i1 = 1
446) i2 = 1
447) endif
448) if (j1 > 0) then
449) j2 = max(min(j1+1,sizej),1)
450) else
451) j1 = 1
452) j2 = 1
453) endif
454) if (i2 == i1) then
455) if (j2 == j1) then
456) ! corner of domain
457) result = this%data(i1+(j1-1)*sizei)
458) else
459) y1 = this%axis2%values(j1)
460) y2 = this%axis2%values(j2)
461) z1 = this%data(i1+(j1-1)*sizei)
462) z3 = this%data(i2+(j2-1)*sizei)
463) call Interpolate(y2,y1,lookup2,z3,z1,result)
464) endif
465) else if (j2 == j1) then
466) x1 = this%axis1%values(i1)
467) x2 = this%axis1%values(i2)
468) z1 = this%data(i1+(j1-1)*sizei)
469) z2 = this%data(i2+(j2-1)*sizei)
470) call Interpolate(x2,x1,lookup1,z2,z1,result)
471) else
472) x1 = this%axis1%values(i1)
473) x2 = this%axis1%values(i2)
474) y1 = this%axis2%values(j1)
475) y2 = this%axis2%values(j2)
476) z1 = this%data(i1+(j1-1)*sizei)
477) z2 = this%data(i2+(j1-1)*sizei)
478) z3 = this%data(i1+(j2-1)*sizei)
479) z4 = this%data(i2+(j2-1)*sizei)
480) result = InterpolateBilinear(lookup1,lookup2,x1,x2,y1,y2,z1,z2,z3,z4)
481) endif
482)
483) end subroutine LookupTableInterpolate2DUniform
484)
485) ! ************************************************************************** !
486)
487) subroutine LookupTableInterpolate2DGeneral(this,lookup1,lookup2,result)
488) !
489) ! Author: Glenn Hammond
490) ! Date: 10/15/14
491) !
492) use Utility_module, only : Interpolate
493)
494) implicit none
495)
496) class(lookup_table_general_type) :: this
497) PetscInt :: i1a
498) PetscInt :: i1b
499) PetscInt :: ja
500) PetscReal :: lookup1
501) PetscReal :: lookup2
502) PetscReal :: result
503)
504) PetscInt :: i2a, i2b
505) PetscInt :: jb
506) PetscInt :: ii1, ii2, iia, iib
507) PetscReal :: x1, x2, z1, z2, xa, xb
508) PetscReal :: interp_a, interp_b
509) PetscInt :: sizei, sizej
510)
511)
512) ! x2,y2,z4
513) ! /|
514) ! / |
515) ! / |
516) ! / |
517) ! / |
518) ! / |
519) ! x1,y2,z3 |
520) ! | |
521) ! | x,y |
522) ! | |
523) ! x1,y1,z1 -x2,y1,z2
524)
525) result = UNINITIALIZED_DOUBLE
526) sizei = this%dims(2)
527) sizej = this%dims(1)
528) ! index axes
529) ja = this%axis1%saved_index
530) i1a = this%axis2%saved_index
531) i1b = this%axis2%saved_index2
532) if (ja > 0) then
533) jb = max(min(ja+1,sizej),1)
534) else
535) ja = 1
536) jb = 1
537) endif
538) if (i1a > 0) then
539) i2a = max(min(i1a+1,sizei),1)
540) else
541) i1a = 1
542) i2a = 1
543) endif
544) if (i1b > 0) then
545) i2b = max(min(i1b+1,sizei),1)
546) else
547) i1b = 1
548) i2b = 1
549) endif
550) if (jb == ja) then
551) ! only use ja/i1a/i2a
552) if (i2a == i1a) then
553) ! corner of domain
554) result = this%data(i1a+(ja-1)*sizei)
555) else
556) ii1 = i1a+(ja-1)*sizei
557) ii2 = i2a+(ja-1)*sizei
558) x1 = this%axis2%values(ii1)
559) x2 = this%axis2%values(ii2)
560) z1 = this%data(ii1)
561) z2 = this%data(ii2)
562) call Interpolate(x2,x1,lookup2,z2,z1,result)
563) endif
564) else
565) ! ja / i*a
566) if (i2a == i1a) then
567) interp_a = this%data(i1a+(ja-1)*sizei)
568) else
569) iia = i1a+(ja-1)*sizei
570) iib = i2a+(ja-1)*sizei
571) x1 = this%axis2%values(iia)
572) x2 = this%axis2%values(iib)
573) z1 = this%data(iia)
574) z2 = this%data(iib)
575) call Interpolate(x2,x1,lookup2,z2,z1,interp_a)
576) endif
577) ! jb / i*b
578) if (i2b == i1b) then
579) interp_b = this%data(i1b+(jb-1)*sizei)
580) else
581) iia = i1b+(jb-1)*sizei
582) iib = i2b+(jb-1)*sizei
583) x1 = this%axis2%values(iia)
584) x2 = this%axis2%values(iib)
585) z1 = this%data(iia)
586) z2 = this%data(iib)
587) call Interpolate(x2,x1,lookup2,z2,z1,interp_b)
588) endif
589) xa = this%axis1%values(ja)
590) xb = this%axis1%values(jb)
591) call Interpolate(xb,xa,lookup1,interp_b,interp_a,result)
592) endif
593)
594) end subroutine LookupTableInterpolate2DGeneral
595)
596) ! ************************************************************************** !
597)
598) subroutine LookupTableTest1D(lookup_table,lookup,desired_result)
599) !
600) ! Deallocates any allocated pointers in axis
601) !
602) ! Author: Glenn Hammond
603) ! Date: 10/15/14
604) !
605) use Utility_module
606)
607) implicit none
608)
609) class(lookup_table_base_type) :: lookup_table
610) PetscReal :: lookup
611) PetscReal :: desired_result
612)
613) PetscReal :: result
614)
615) 100 format(3(f12.6),l)
616)
617) result = lookup_table%Sample(lookup)
618) write(*,100) lookup, desired_result, result, Equal(result,desired_result)
619)
620) end subroutine LookupTableTest1D
621)
622) ! ************************************************************************** !
623)
624) subroutine LookupTableTest2D(lookup_table,lookup1,lookup2,desired_result)
625) !
626) ! Deallocates any allocated pointers in axis
627) !
628) ! Author: Glenn Hammond
629) ! Date: 10/15/14
630) !
631) use Utility_module
632)
633) implicit none
634)
635) class(lookup_table_base_type) :: lookup_table
636) PetscReal :: lookup1
637) PetscReal :: lookup2
638) PetscReal :: desired_result
639)
640) PetscReal :: result
641)
642) 100 format(4(f12.6),l)
643)
644) result = lookup_table%Sample(lookup1,lookup2)
645) write(*,100) lookup1, lookup2, desired_result, result, &
646) Equal(result,desired_result)
647)
648) end subroutine LookupTableTest2D
649)
650) ! ************************************************************************** !
651)
652) subroutine LookupTableAxisDestroy(axis)
653) !
654) ! Deallocates any allocated pointers in axis
655) !
656) ! Author: Glenn Hammond
657) ! Date: 10/15/14
658) !
659) use Utility_module
660)
661) implicit none
662)
663) class(lookup_table_axis_type), pointer :: axis
664)
665) if (.not.associated(axis)) return
666)
667) call DeallocateArray(axis%values)
668) deallocate(axis)
669) nullify(axis)
670)
671) end subroutine LookupTableAxisDestroy
672)
673) ! ************************************************************************** !
674)
675) subroutine LookupTableUniformDestroy(lookup_table)
676) !
677) ! Deallocates any allocated pointers in lookup table
678) !
679) ! Author: Glenn Hammond
680) ! Date: 10/15/14
681) !
682) use Utility_module
683)
684) implicit none
685)
686) class(lookup_table_uniform_type), pointer :: lookup_table
687)
688) call LookupTableAxisDestroy(lookup_table%axis1)
689) call LookupTableAxisDestroy(lookup_table%axis2)
690) call LookupTableAxisDestroy(lookup_table%axis3)
691) call DeallocateArray(lookup_table%data)
692) deallocate(lookup_table)
693) nullify(lookup_table)
694)
695) end subroutine LookupTableUniformDestroy
696)
697) ! ************************************************************************** !
698)
699) subroutine LookupTableGeneralDestroy(lookup_table)
700) !
701) ! Deallocates any allocated pointers in lookup table
702) !
703) ! Author: Glenn Hammond
704) ! Date: 10/15/14
705) !
706) use Utility_module
707)
708) implicit none
709)
710) class(lookup_table_general_type), pointer :: lookup_table
711)
712) if (.not.associated(lookup_table)) return
713)
714) call LookupTableAxisDestroy(lookup_table%axis1)
715) ! axis2 is different type
716) if (associated(lookup_table%axis2)) then
717) call DeallocateArray(lookup_table%axis2%values)
718) deallocate(lookup_table%axis2)
719) nullify(lookup_table%axis2)
720) endif
721) call DeallocateArray(lookup_table%data)
722) deallocate(lookup_table)
723) nullify(lookup_table)
724)
725) end subroutine LookupTableGeneralDestroy
726)
727) end module Lookup_Table_module