V3FIT
equilibrium.f
Go to the documentation of this file.
1 !*******************************************************************************
4 !
5 ! Note separating the Doxygen comment block here so detailed decription is
6 ! found in the Module not the file.
7 !
17 !*******************************************************************************
18 
19  MODULE equilibrium
22 
23  IMPLICIT NONE
24 
25 !*******************************************************************************
26 ! equilibrium module parameters
27 !*******************************************************************************
29  INTEGER, PARAMETER :: equilibrium_no_type = -1
32  INTEGER, PARAMETER :: equilibrium_vmec_type = 0
35  INTEGER, PARAMETER :: equilibrium_vacuum_type = 1
38  INTEGER, PARAMETER :: equilibrium_siesta_type = 2
39 
40 !*******************************************************************************
41 ! DERIVED-TYPE DECLARATIONS
42 ! 1) equilibrium base class
43 !
44 !*******************************************************************************
45 !-------------------------------------------------------------------------------
49 !-------------------------------------------------------------------------------
55  INTEGER :: type = equilibrium_no_type
57  LOGICAL :: force_solve
59  TYPE (vmec_class), POINTER :: vmec => null()
61  TYPE (vacuum_class), POINTER :: vacuum => null()
63  TYPE (siesta_class), POINTER :: siesta => null()
64  END TYPE
65 
66 !*******************************************************************************
67 ! INTERFACE BLOCKS
68 !*******************************************************************************
69 !-------------------------------------------------------------------------------
73 !-------------------------------------------------------------------------------
75  MODULE PROCEDURE equilibrium_construct_vmec, &
78  END INTERFACE
79 
80 !-------------------------------------------------------------------------------
82 !-------------------------------------------------------------------------------
84  MODULE PROCEDURE equilibrium_set_magnetic_cache_response, &
86  END INTERFACE
87 
88 !-------------------------------------------------------------------------------
90 !-------------------------------------------------------------------------------
92  MODULE PROCEDURE equilibrium_get_ne_cart, &
94  END INTERFACE
95 
96 !-------------------------------------------------------------------------------
98 !-------------------------------------------------------------------------------
100  MODULE PROCEDURE equilibrium_get_gp_ne_ij, &
103  END INTERFACE
104 
105 !-------------------------------------------------------------------------------
107 !-------------------------------------------------------------------------------
109  MODULE PROCEDURE equilibrium_get_te_cart, &
111  END INTERFACE
112 
113 !-------------------------------------------------------------------------------
116 !-------------------------------------------------------------------------------
118  MODULE PROCEDURE equilibrium_get_gp_te_ij, &
121  END INTERFACE
122 
123 !-------------------------------------------------------------------------------
125 !-------------------------------------------------------------------------------
127  MODULE PROCEDURE equilibrium_get_ti_cart, &
129  END INTERFACE
130 
131 !-------------------------------------------------------------------------------
134 !-------------------------------------------------------------------------------
136  MODULE PROCEDURE equilibrium_get_gp_ti_ij, &
139  END INTERFACE
140 
141 !-------------------------------------------------------------------------------
143 !-------------------------------------------------------------------------------
145  MODULE PROCEDURE equilibrium_get_ze_cart, &
147  END INTERFACE
148 
149 !-------------------------------------------------------------------------------
151 !-------------------------------------------------------------------------------
153  MODULE PROCEDURE equilibrium_get_sxrem_cart, &
155  END INTERFACE
156 
157 !-------------------------------------------------------------------------------
160 !-------------------------------------------------------------------------------
162  MODULE PROCEDURE equilibrium_get_gp_sxrem_ij, &
165  END INTERFACE
166 
167 !-------------------------------------------------------------------------------
169 !-------------------------------------------------------------------------------
171  MODULE PROCEDURE equilibrium_get_p_cart, &
173  END INTERFACE
174 
175  CONTAINS
176 !*******************************************************************************
177 ! CONSTRUCTION SUBROUTINES
178 !*******************************************************************************
179 !-------------------------------------------------------------------------------
188 !-------------------------------------------------------------------------------
189  FUNCTION equilibrium_construct_vmec(vmec_object, force_solve)
190 
191  IMPLICIT NONE
192 
193 ! Declare Arguments
195  TYPE (vmec_class), POINTER :: vmec_object
196  LOGICAL, INTENT(in) :: force_solve
197 
198 ! local variables
199  REAL (rprec) :: start_time
200 
201 ! Start of executable code
202  start_time = profiler_get_start_time()
203 
205 
207  equilibrium_construct_vmec%vmec => vmec_object
208  equilibrium_construct_vmec%force_solve = force_solve
209 
210  CALL profiler_set_stop_time('equilibrium_construct_vmec', &
211  & start_time)
212 
213  END FUNCTION
214 
215 !-------------------------------------------------------------------------------
225 !-------------------------------------------------------------------------------
226  FUNCTION equilibrium_construct_vacuum(vacuum_object, force_solve)
227 
228  IMPLICIT NONE
229 
230 ! Declare Arguments
232  TYPE (vacuum_class), POINTER :: vacuum_object
233  LOGICAL, INTENT(in) :: force_solve
234 
235 ! local variables
236  REAL (rprec) :: start_time
237 
238 ! Start of executable code
239  start_time = profiler_get_start_time()
240 
242 
244  equilibrium_construct_vacuum%vacuum => vacuum_object
245  equilibrium_construct_vacuum%force_solve = force_solve
246 
247  CALL profiler_set_stop_time('equilibrium_construct_vacuum', &
248  & start_time)
249 
250  END FUNCTION
251 
252 !-------------------------------------------------------------------------------
262 !-------------------------------------------------------------------------------
263  FUNCTION equilibrium_construct_siesta(siesta_object, force_solve)
264 
265  IMPLICIT NONE
266 
267 ! Declare Arguments
269  TYPE (siesta_class), POINTER :: siesta_object
270  LOGICAL, INTENT(in) :: force_solve
271 
272 ! local variables
273  REAL (rprec) :: start_time
274 
275 ! Start of executable code
276  start_time = profiler_get_start_time()
277 
279 
281  equilibrium_construct_siesta%siesta => siesta_object
282  equilibrium_construct_siesta%force_solve = force_solve
283 
284  CALL profiler_set_stop_time('equilibrium_construct_siesta', &
285  & start_time)
286 
287  END FUNCTION
288 
289 !*******************************************************************************
290 ! DESTRUCTION SUBROUTINES
291 !*******************************************************************************
292 !-------------------------------------------------------------------------------
298 !-------------------------------------------------------------------------------
299  SUBROUTINE equilibrium_destruct(this)
300 
301  IMPLICIT NONE
302 
303 ! Declare Arguments
304  TYPE (equilibrium_class), POINTER :: this
305 
306 ! Start of executable code
307  IF (ASSOCIATED(this%vmec)) THEN
308  CALL vmec_destruct(this%vmec)
309  this%vmec => null()
310  END IF
311 
312  IF (ASSOCIATED(this%vacuum)) THEN
313  CALL vacuum_destruct(this%vacuum)
314  this%vacuum => null()
315  END IF
316 
317  IF (ASSOCIATED(this%siesta)) THEN
318  CALL siesta_destruct(this%siesta)
319  this%siesta => null()
320  END IF
321 
322  this%type = equilibrium_no_type
323  this%force_solve = .false.
324 
325  DEALLOCATE(this)
326 
327  END SUBROUTINE
328 
329 !*******************************************************************************
330 ! SETTER SUBROUTINES
331 !*******************************************************************************
332 !-------------------------------------------------------------------------------
350 !-------------------------------------------------------------------------------
351  SUBROUTINE equilibrium_set_param(this, id, i_index, j_index, &
352  & value, eq_comm, state_flags)
353 
354  IMPLICIT NONE
355 
356 ! Declare Arguments
357  TYPE (equilibrium_class), INTENT(inout) :: this
358  INTEGER, INTENT(in) :: id
359  INTEGER, INTENT(in) :: i_index
360  INTEGER, INTENT(in) :: j_index
361  REAL (rprec), INTENT(in) :: value
362  INTEGER, INTENT(in) :: eq_comm
363  INTEGER, INTENT(inout) :: state_flags
364 
365 ! local variables
366  REAL (rprec) :: start_time
367 
368 ! Start of executable code
369  start_time = profiler_get_start_time()
370 
371  SELECT CASE (this%type)
372 
373 ! Set param needs to come first so that function gets executed. Otherwise, when
374 ! this%not_converged was true, the expression will evaluate to false all the
375 ! time and not bother calling vmec_set_param.
376  CASE (equilibrium_vmec_type)
377  CALL vmec_set_param(this%vmec, id, i_index, j_index, value, &
378  & eq_comm, state_flags)
379 
381  CALL vacuum_set_param(this%vacuum, id, i_index, value)
382 
384  CALL siesta_set_param(this%siesta, id, i_index, j_index, &
385  & value, eq_comm, state_flags)
386 
387  END SELECT
388 
389  CALL profiler_set_stop_time('equilibrium_set_param', start_time)
390 
391  END SUBROUTINE
392 
393 !-------------------------------------------------------------------------------
410 !-------------------------------------------------------------------------------
411  SUBROUTINE equilibrium_set_magnetic_cache_response(this, &
412  & response_object, &
413  & state_flags)
415  USE model_state
416 
417  IMPLICIT NONE
418 
419 ! Declare Arguments
420  TYPE (equilibrium_class), INTENT(inout) :: this
421  TYPE (magnetic_response_class), INTENT(in) :: response_object
422  INTEGER, INTENT(in) :: state_flags
423 
424 ! local variables
425  REAL (rprec) :: start_time
426 
427 ! Start of executable code
428  start_time = profiler_get_start_time()
429 
430  SELECT CASE (this%type)
431 
432  CASE (equilibrium_vmec_type)
433  CALL vmec_set_magnetic_cache(this%vmec, response_object)
434 
435 ! If the equilibrium is already converged, compute the magnetic cache as well.
436  IF (.not.btest(state_flags, model_state_vmec_flag)) THEN
437  CALL vmec_set_magnetic_cache(this%vmec)
438  END IF
439 
441  CALL siesta_set_magnetic_cache(this%siesta, response_object)
442 
443 ! If the equilibrium is already converged, compute the magnetic cache as well.
444  IF (.not.btest(state_flags, model_state_siesta_flag)) THEN
445  CALL siesta_set_magnetic_cache(this%siesta)
446  END IF
447 
448  END SELECT
449 
450  CALL profiler_set_stop_time( &
451  & 'equilibrium_set_magnetic_cache_response', start_time)
452 
453  END SUBROUTINE
454 
455 !-------------------------------------------------------------------------------
472 !-------------------------------------------------------------------------------
473  SUBROUTINE equilibrium_set_magnetic_cache_point(this, use_axi, &
474  & state_flags)
476 
477  IMPLICIT NONE
478 
479 ! Declare Arguments
480  TYPE (equilibrium_class), INTENT(inout) :: this
481  LOGICAL, INTENT(in) :: use_axi
482  INTEGER, INTENT(in) :: state_flags
483 
484 ! local variables
485  REAL (rprec) :: start_time
486 
487 ! Start of executable code
488  start_time = profiler_get_start_time()
489 
490  SELECT CASE (this%type)
491 
492  CASE (equilibrium_vmec_type)
493  CALL vmec_set_magnetic_cache(this%vmec, use_axi)
494 
495 ! If the equilibrium is already converged, compute the magnetic cache as well.
496  IF (.not.btest(state_flags, model_state_vmec_flag)) THEN
497  CALL vmec_set_magnetic_cache(this%vmec)
498  END IF
499 
501  CALL siesta_set_magnetic_cache(this%siesta, use_axi)
502 
503 ! If the equilibrium is already converged, compute the magnetic cache as well.
504  IF (.not.btest(state_flags, model_state_siesta_flag)) THEN
505  CALL siesta_set_magnetic_cache(this%siesta)
506  END IF
507 
508  END SELECT
509 
510  CALL profiler_set_stop_time( &
511  & 'equilibrium_set_magnetic_cache_point', start_time)
512 
513  END SUBROUTINE
514 
515 !*******************************************************************************
516 ! GETTER SUBROUTINES
517 !*******************************************************************************
518 !-------------------------------------------------------------------------------
530 !-------------------------------------------------------------------------------
531  FUNCTION equilibrium_get_param_id(this, param_name)
532 
533  IMPLICIT NONE
534 
535 ! Declare Arguments
536  INTEGER :: equilibrium_get_param_id
537  TYPE (equilibrium_class), INTENT(in) :: this
538  CHARACTER (len=*), INTENT(in) :: param_name
539 
540 ! local variables
541  REAL (rprec) :: start_time
542 
543 ! Start of executable code
544  start_time = profiler_get_start_time()
545 
546  SELECT CASE (this%type)
547 
548  CASE (equilibrium_vmec_type)
550  & vmec_get_param_id(this%vmec, param_name)
551 
554  & vacuum_get_param_id(this%vacuum, param_name)
555 
558  & siesta_get_param_id(this%siesta, param_name)
559 
560  CASE DEFAULT
561  equilibrium_get_param_id = data_no_id
562 
563  END SELECT
564 
565  CALL profiler_set_stop_time('equilibrium_get_param_id', &
566  & start_time)
567 
568  END FUNCTION
569 
570 !-------------------------------------------------------------------------------
584 !-------------------------------------------------------------------------------
585  FUNCTION equilibrium_get_param_value(this, id, i_index, j_index)
586 
587  IMPLICIT NONE
588 
589 ! Declare Arguments
590  REAL (rprec) :: equilibrium_get_param_value
591  TYPE (equilibrium_class), INTENT(in) :: this
592  INTEGER, INTENT(in) :: id
593  INTEGER, INTENT(in) :: i_index
594  INTEGER, INTENT(in) :: j_index
595 
596 ! local variables
597  REAL (rprec) :: start_time
598 
599 ! Start of executable code
600  start_time = profiler_get_start_time()
601 
602  SELECT CASE (this%type)
603 
604  CASE (equilibrium_vmec_type)
606  & vmec_get_param_value(this%vmec, id, i_index, j_index)
607 
610  & vacuum_get_param_value(this%vacuum, id, i_index)
611 
614  & siesta_get_param_value(this%siesta, id, i_index, j_index)
615 
616  CASE DEFAULT
618 
619  END SELECT
620 
621  CALL profiler_set_stop_time('equilibrium_get_param_value', &
622  & start_time)
623 
624  END FUNCTION
625 
626 !-------------------------------------------------------------------------------
638 !-------------------------------------------------------------------------------
639  FUNCTION equilibrium_get_param_name(this, id)
640 
641  IMPLICIT NONE
642 
643 ! Declare Arguments
644  CHARACTER(len=data_name_length) :: equilibrium_get_param_name
645  TYPE (equilibrium_class), INTENT(in) :: this
646  INTEGER, INTENT(in) :: id
647 
648 ! local variables
649  REAL (rprec) :: start_time
650 
651 ! Start of executable code
652  start_time = profiler_get_start_time()
653 
654  SELECT CASE (this%type)
655 
656  CASE (equilibrium_vmec_type)
658  & vmec_get_param_name(this%vmec, id)
659 
662  & vacuum_get_param_name(this%vacuum, id)
663 
666  & siesta_get_param_name(this%siesta, id)
667 
668  CASE DEFAULT
670 
671  END SELECT
672 
673  CALL profiler_set_stop_time('equilibrium_get_param_name', &
674  & start_time)
675 
676  END FUNCTION
677 
678 !-------------------------------------------------------------------------------
688 !-------------------------------------------------------------------------------
690 
691  IMPLICIT NONE
692 
693 ! Declare Arguments
695  TYPE (equilibrium_class), INTENT(in) :: this
696 
697 ! local variables
698  REAL (rprec) :: start_time
699 
700 ! Start of executable code
701  start_time = profiler_get_start_time()
702 
703  SELECT CASE (this%type)
704 
705  CASE (equilibrium_vmec_type)
707  & vmec_get_gp_ne_num_hyper_param(this%vmec)
708 
711  & siesta_get_gp_ne_num_hyper_param(this%siesta)
712 
713  CASE DEFAULT
715 
716  END SELECT
717 
718  CALL profiler_set_stop_time( &
719  & 'equilibrium_get_gp_ne_num_hyper_param', start_time)
720 
721  END FUNCTION
722 
723 !-------------------------------------------------------------------------------
733 !-------------------------------------------------------------------------------
734  FUNCTION equilibrium_get_ne_af(this)
735 
736  IMPLICIT NONE
737 
738 ! Declare Arguments
739  REAL (rprec), DIMENSION(:), POINTER :: equilibrium_get_ne_af
740  TYPE (equilibrium_class), INTENT(in) :: this
741 
742 ! local variables
743  REAL (rprec) :: start_time
744 
745 ! Start of executable code
746  start_time = profiler_get_start_time()
747 
748  SELECT CASE (this%type)
749 
750  CASE (equilibrium_vmec_type)
751  equilibrium_get_ne_af => vmec_get_ne_af(this%vmec)
752 
755 
756  CASE DEFAULT
757  equilibrium_get_ne_af => null()
758 
759  END SELECT
760 
761  END FUNCTION
762 
763 !-------------------------------------------------------------------------------
777 !-------------------------------------------------------------------------------
778  FUNCTION equilibrium_get_gp_ne_ij(this, i, j)
779 
780  IMPLICIT NONE
781 
782 ! Declare Arguments
783  REAL (rprec) :: equilibrium_get_gp_ne_ij
784  TYPE (equilibrium_class), INTENT(in) :: this
785  INTEGER, INTENT(in) :: i
786  INTEGER, INTENT(in) :: j
787 
788 ! local variables
789  REAL (rprec) :: start_time
790 
791 ! Start of executable code
792  start_time = profiler_get_start_time()
793 
794  SELECT CASE (this%type)
795 
796  CASE (equilibrium_vmec_type)
797  equilibrium_get_gp_ne_ij = vmec_get_gp_ne(this%vmec, i, j)
798 
801  & i, j)
802 
803  CASE DEFAULT
805 
806  END SELECT
807 
808  CALL profiler_set_stop_time('equilibrium_get_gp_ne_ij', &
809  & start_time)
810 
811  END FUNCTION
812 
813 !-------------------------------------------------------------------------------
827 !-------------------------------------------------------------------------------
828  FUNCTION equilibrium_get_gp_ne_pi(this, x_cart, i)
829 
830  IMPLICIT NONE
831 
832 ! Declare Arguments
833  REAL (rprec) :: equilibrium_get_gp_ne_pi
834  TYPE (equilibrium_class), INTENT(in) :: this
835  REAL (rprec), DIMENSION(3), INTENT(in) :: x_cart
836  INTEGER, INTENT(in) :: i
837 
838 ! local variables
839  REAL (rprec) :: start_time
840 
841 ! Start of executable code
842  start_time = profiler_get_start_time()
843 
844  SELECT CASE (this%type)
845 
846  CASE (equilibrium_vmec_type)
847  equilibrium_get_gp_ne_pi = vmec_get_gp_ne(this%vmec, &
848  & x_cart, i)
849 
852  & x_cart, i)
853 
854  CASE DEFAULT
856 
857  END SELECT
858 
859  CALL profiler_set_stop_time('equilibrium_get_gp_ne_pi', &
860  & start_time)
861 
862  END FUNCTION
863 
864 !-------------------------------------------------------------------------------
879 !-------------------------------------------------------------------------------
880  FUNCTION equilibrium_get_gp_ne_pp(this, x_cart, y_cart)
881 
882  IMPLICIT NONE
883 
884 ! Declare Arguments
885  REAL (rprec) :: equilibrium_get_gp_ne_pp
886  TYPE (equilibrium_class), INTENT(in) :: this
887  REAL (rprec), DIMENSION(3), INTENT(in) :: x_cart
888  REAL (rprec), DIMENSION(3), INTENT(in) :: y_cart
889 
890 ! local variables
891  REAL (rprec) :: start_time
892 
893 ! Start of executable code
894  start_time = profiler_get_start_time()
895 
896  SELECT CASE (this%type)
897 
898  CASE (equilibrium_vmec_type)
899  equilibrium_get_gp_ne_pp = vmec_get_gp_ne(this%vmec, &
900  & x_cart, y_cart)
901 
904  & x_cart, y_cart)
905 
906  CASE DEFAULT
908 
909  END SELECT
910 
911  CALL profiler_set_stop_time('equilibrium_get_gp_ne_pp', &
912  & start_time)
913 
914  END FUNCTION
915 
916 !-------------------------------------------------------------------------------
929 !-------------------------------------------------------------------------------
930  FUNCTION equilibrium_get_ne_cart(this, x_cart)
931 
932  IMPLICIT NONE
933 
934 ! Declare Arguments
935  REAL (rprec) :: equilibrium_get_ne_cart
936  TYPE (equilibrium_class), INTENT(in) :: this
937  REAL (rprec), DIMENSION(3), INTENT(in) :: x_cart
938 
939 ! local variables
940  REAL (rprec) :: start_time
941 
942 ! Start of executable code
943  start_time = profiler_get_start_time()
944 
945  SELECT CASE (this%type)
946 
947  CASE (equilibrium_vmec_type)
948  equilibrium_get_ne_cart = vmec_get_ne(this%vmec, x_cart)
949 
951  equilibrium_get_ne_cart = siesta_get_ne(this%siesta, x_cart)
952 
953  CASE DEFAULT
955 
956  END SELECT
957 
958  CALL profiler_set_stop_time('equilibrium_get_ne_cart', &
959  & start_time)
960 
961  END FUNCTION
962 
963 !-------------------------------------------------------------------------------
976 !-------------------------------------------------------------------------------
977  FUNCTION equilibrium_get_ne_radial(this, r)
978 
979  IMPLICIT NONE
980 
981 ! Declare Arguments
982  REAL (rprec) :: equilibrium_get_ne_radial
983  TYPE (equilibrium_class), INTENT(in) :: this
984  REAL (rprec), INTENT(in) :: r
985 
986 ! local variables
987  REAL (rprec) :: start_time
988 
989 ! Start of executable code
990  start_time = profiler_get_start_time()
991 
992  SELECT CASE (this%type)
993 
994  CASE (equilibrium_vmec_type)
995  equilibrium_get_ne_radial = vmec_get_ne(this%vmec, r)
996 
998  equilibrium_get_ne_radial = siesta_get_ne(this%siesta, r)
999 
1000  CASE DEFAULT
1002 
1003  END SELECT
1004 
1005  CALL profiler_set_stop_time('equilibrium_get_ne_radial', &
1006  & start_time)
1007 
1008  END FUNCTION
1009 
1010 !-------------------------------------------------------------------------------
1020 !-------------------------------------------------------------------------------
1023  IMPLICIT NONE
1024 
1025 ! Declare Arguments
1027  TYPE (equilibrium_class), INTENT(in) :: this
1028 
1029 ! local variables
1030  REAL (rprec) :: start_time
1031 
1032 ! Start of executable code
1033  start_time = profiler_get_start_time()
1034 
1035  SELECT CASE (this%type)
1036 
1037  CASE (equilibrium_vmec_type)
1039  & vmec_get_gp_te_num_hyper_param(this%vmec)
1040 
1043  & siesta_get_gp_te_num_hyper_param(this%siesta)
1044 
1045  CASE DEFAULT
1047 
1048  END SELECT
1049 
1050  CALL profiler_set_stop_time( &
1051  & 'equilibrium_get_gp_te_num_hyper_param', start_time)
1052 
1053  END FUNCTION
1054 
1055 !-------------------------------------------------------------------------------
1065 !-------------------------------------------------------------------------------
1066  FUNCTION equilibrium_get_te_af(this)
1068  IMPLICIT NONE
1069 
1070 ! Declare Arguments
1071  REAL (rprec), DIMENSION(:), POINTER :: equilibrium_get_te_af
1072  TYPE (equilibrium_class), INTENT(in) :: this
1073 
1074 ! local variables
1075  REAL (rprec) :: start_time
1076 
1077 ! Start of executable code
1078  start_time = profiler_get_start_time()
1079 
1080  SELECT CASE (this%type)
1081 
1082  CASE (equilibrium_vmec_type)
1083  equilibrium_get_te_af => vmec_get_te_af(this%vmec)
1084 
1086  equilibrium_get_te_af => siesta_get_te_af(this%siesta)
1087 
1088  CASE DEFAULT
1089  equilibrium_get_te_af => null()
1090 
1091  END SELECT
1092 
1093  CALL profiler_set_stop_time('equilibrium_get_te_af', start_time)
1094 
1095  END FUNCTION
1096 
1097 !-------------------------------------------------------------------------------
1111 !-------------------------------------------------------------------------------
1112  FUNCTION equilibrium_get_gp_te_ij(this, i, j)
1114  IMPLICIT NONE
1115 
1116 ! Declare Arguments
1117  REAL (rprec) :: equilibrium_get_gp_te_ij
1118  TYPE (equilibrium_class), INTENT(in) :: this
1119  INTEGER, INTENT(in) :: i
1120  INTEGER, INTENT(in) :: j
1121 
1122 ! local variables
1123  REAL (rprec) :: start_time
1124 
1125 ! Start of executable code
1126  start_time = profiler_get_start_time()
1127 
1128  SELECT CASE (this%type)
1129 
1130  CASE (equilibrium_vmec_type)
1131  equilibrium_get_gp_te_ij = vmec_get_gp_te(this%vmec, i, j)
1132 
1135  & i, j)
1136 
1137  CASE DEFAULT
1139 
1140  END SELECT
1141 
1142  CALL profiler_set_stop_time('equilibrium_get_gp_te_ij', &
1143  & start_time)
1144 
1145  END FUNCTION
1146 
1147 !-------------------------------------------------------------------------------
1162 !-------------------------------------------------------------------------------
1163  FUNCTION equilibrium_get_gp_te_pi(this, x_cart, i)
1165  IMPLICIT NONE
1166 
1167 ! Declare Arguments
1168  REAL (rprec) :: equilibrium_get_gp_te_pi
1169  TYPE (equilibrium_class), INTENT(in) :: this
1170  REAL (rprec), DIMENSION(3), INTENT(in) :: x_cart
1171  INTEGER, INTENT(in) :: i
1172 
1173 ! local variables
1174  REAL (rprec) :: start_time
1175 
1176 ! Start of executable code
1177  start_time = profiler_get_start_time()
1178 
1179  SELECT CASE (this%type)
1180 
1181  CASE (equilibrium_vmec_type)
1182  equilibrium_get_gp_te_pi = vmec_get_gp_te(this%vmec, &
1183  & x_cart, i)
1184 
1187  & x_cart, i)
1188 
1189  CASE DEFAULT
1191 
1192  END SELECT
1193 
1194  CALL profiler_set_stop_time('equilibrium_get_gp_te_pi', &
1195  & start_time)
1196 
1197  END FUNCTION
1198 
1199 !-------------------------------------------------------------------------------
1214 !-------------------------------------------------------------------------------
1215  FUNCTION equilibrium_get_gp_te_pp(this, x_cart, y_cart)
1217  IMPLICIT NONE
1218 
1219 ! Declare Arguments
1220  REAL (rprec) :: equilibrium_get_gp_te_pp
1221  TYPE (equilibrium_class), INTENT(in) :: this
1222  REAL (rprec), DIMENSION(3), INTENT(in) :: x_cart
1223  REAL (rprec), DIMENSION(3), INTENT(in) :: y_cart
1224 
1225 ! local variables
1226  REAL (rprec) :: start_time
1227 
1228 ! Start of executable code
1229  start_time = profiler_get_start_time()
1230 
1231  SELECT CASE (this%type)
1232 
1233  CASE (equilibrium_vmec_type)
1234  equilibrium_get_gp_te_pp = vmec_get_gp_te(this%vmec, &
1235  & x_cart, y_cart)
1236 
1239  & x_cart, y_cart)
1240 
1241  CASE DEFAULT
1243 
1244  END SELECT
1245 
1246  CALL profiler_set_stop_time('equilibrium_get_gp_te_pp', &
1247  & start_time)
1248 
1249  END FUNCTION
1250 
1251 !-------------------------------------------------------------------------------
1264 !-------------------------------------------------------------------------------
1265  FUNCTION equilibrium_get_te_cart(this, x_cart)
1267  IMPLICIT NONE
1268 
1269 ! Declare Arguments
1270  REAL (rprec) :: equilibrium_get_te_cart
1271  TYPE (equilibrium_class), INTENT(in) :: this
1272  REAL (rprec), DIMENSION(3), INTENT(in) :: x_cart
1273 
1274 ! local variables
1275  REAL (rprec) :: start_time
1276 
1277 ! Start of executable code
1278  start_time = profiler_get_start_time()
1279 
1280  SELECT CASE (this%type)
1281 
1282  CASE (equilibrium_vmec_type)
1283  equilibrium_get_te_cart = vmec_get_te(this%vmec, x_cart)
1284 
1286  equilibrium_get_te_cart = siesta_get_te(this%siesta, x_cart)
1287 
1288  CASE DEFAULT
1290 
1291  END SELECT
1292 
1293  CALL profiler_set_stop_time('equilibrium_get_te_cart', &
1294  & start_time)
1295 
1296  END FUNCTION
1297 
1298 !-------------------------------------------------------------------------------
1311 !-------------------------------------------------------------------------------
1312  FUNCTION equilibrium_get_te_radial(this, r)
1314  IMPLICIT NONE
1315 
1316 ! Declare Arguments
1317  REAL (rprec) :: equilibrium_get_te_radial
1318  TYPE (equilibrium_class), INTENT(in) :: this
1319  REAL (rprec), INTENT(in) :: r
1320 
1321 ! local variables
1322  REAL (rprec) :: start_time
1323 
1324 ! Start of executable code
1325  start_time = profiler_get_start_time()
1326 
1327  SELECT CASE (this%type)
1328 
1329  CASE (equilibrium_vmec_type)
1330  equilibrium_get_te_radial = vmec_get_te(this%vmec, r)
1331 
1333  equilibrium_get_te_radial = siesta_get_te(this%siesta, r)
1334 
1335  CASE DEFAULT
1337 
1338  END SELECT
1339 
1340  CALL profiler_set_stop_time('equilibrium_get_te_radial', &
1341  & start_time)
1342 
1343  END FUNCTION
1344 
1345 !-------------------------------------------------------------------------------
1355 !-------------------------------------------------------------------------------
1358  IMPLICIT NONE
1359 
1360 ! Declare Arguments
1362  TYPE (equilibrium_class), INTENT(in) :: this
1363 
1364 ! local variables
1365  REAL (rprec) :: start_time
1366 
1367 ! Start of executable code
1368  start_time = profiler_get_start_time()
1369 
1370  SELECT CASE (this%type)
1371 
1372  CASE (equilibrium_vmec_type)
1374  & vmec_get_gp_ti_num_hyper_param(this%vmec)
1375 
1378  & siesta_get_gp_ti_num_hyper_param(this%siesta)
1379 
1380  CASE DEFAULT
1382 
1383  END SELECT
1384 
1385  CALL profiler_set_stop_time( &
1386  & 'equilibrium_get_gp_ti_num_hyper_param', start_time)
1387 
1388  END FUNCTION
1389 
1390 !-------------------------------------------------------------------------------
1400 !-------------------------------------------------------------------------------
1401  FUNCTION equilibrium_get_ti_af(this)
1403  IMPLICIT NONE
1404 
1405 ! Declare Arguments
1406  REAL (rprec), DIMENSION(:), POINTER :: equilibrium_get_ti_af
1407  TYPE (equilibrium_class), INTENT(in) :: this
1408 
1409 ! local variables
1410  REAL (rprec) :: start_time
1411 
1412 ! Start of executable code
1413  start_time = profiler_get_start_time()
1414 
1415  SELECT CASE (this%type)
1416 
1417  CASE (equilibrium_vmec_type)
1418  equilibrium_get_ti_af => vmec_get_ti_af(this%vmec)
1419 
1421  equilibrium_get_ti_af => siesta_get_ti_af(this%siesta)
1422 
1423  CASE DEFAULT
1424  equilibrium_get_ti_af => null()
1425 
1426  END SELECT
1427 
1428  CALL profiler_set_stop_time('equilibrium_get_ti_af', start_time)
1429 
1430  END FUNCTION
1431 
1432 !-------------------------------------------------------------------------------
1446 !-------------------------------------------------------------------------------
1447  FUNCTION equilibrium_get_gp_ti_ij(this, i, j)
1449  IMPLICIT NONE
1450 
1451 ! Declare Arguments
1452  REAL (rprec) :: equilibrium_get_gp_ti_ij
1453  TYPE (equilibrium_class), INTENT(in) :: this
1454  INTEGER, INTENT(in) :: i
1455  INTEGER, INTENT(in) :: j
1456 
1457 ! local variables
1458  REAL (rprec) :: start_time
1459 
1460 ! Start of executable code
1461  start_time = profiler_get_start_time()
1462 
1463  SELECT CASE (this%type)
1464 
1465  CASE (equilibrium_vmec_type)
1466  equilibrium_get_gp_ti_ij = vmec_get_gp_ti(this%vmec, i, j)
1467 
1470  & i, j)
1471 
1472  CASE DEFAULT
1474 
1475  END SELECT
1476 
1477  CALL profiler_set_stop_time('equilibrium_get_gp_ti_ij', &
1478  & start_time)
1479 
1480  END FUNCTION
1481 
1482 !-------------------------------------------------------------------------------
1496 !-------------------------------------------------------------------------------
1497  FUNCTION equilibrium_get_gp_ti_pi(this, x_cart, i)
1499  IMPLICIT NONE
1500 
1501 ! Declare Arguments
1502  REAL (rprec) :: equilibrium_get_gp_ti_pi
1503  TYPE (equilibrium_class), INTENT(in) :: this
1504  REAL (rprec), DIMENSION(3), INTENT(in) :: x_cart
1505  INTEGER, INTENT(in) :: i
1506 
1507 ! local variables
1508  REAL (rprec) :: start_time
1509 
1510 ! Start of executable code
1511  start_time = profiler_get_start_time()
1512 
1513  SELECT CASE (this%type)
1514 
1515  CASE (equilibrium_vmec_type)
1516  equilibrium_get_gp_ti_pi = vmec_get_gp_ti(this%vmec, &
1517  & x_cart, i)
1518 
1521  & x_cart, i)
1522 
1523  CASE DEFAULT
1525 
1526  END SELECT
1527 
1528  CALL profiler_set_stop_time('equilibrium_get_gp_ti_pi', &
1529  & start_time)
1530 
1531  END FUNCTION
1532 
1533 !-------------------------------------------------------------------------------
1548 !-------------------------------------------------------------------------------
1549  FUNCTION equilibrium_get_gp_ti_pp(this, x_cart, y_cart)
1551  IMPLICIT NONE
1552 
1553 ! Declare Arguments
1554  REAL (rprec) :: equilibrium_get_gp_ti_pp
1555  TYPE (equilibrium_class), INTENT(in) :: this
1556  REAL (rprec), DIMENSION(3), INTENT(in) :: x_cart
1557  REAL (rprec), DIMENSION(3), INTENT(in) :: y_cart
1558 
1559 ! local variables
1560  REAL (rprec) :: start_time
1561 
1562 ! Start of executable code
1563  start_time = profiler_get_start_time()
1564 
1565  SELECT CASE (this%type)
1566 
1567  CASE (equilibrium_vmec_type)
1568  equilibrium_get_gp_ti_pp = vmec_get_gp_ti(this%vmec, &
1569  & x_cart, y_cart)
1570 
1573  & x_cart, y_cart)
1574 
1575  CASE DEFAULT
1577 
1578  END SELECT
1579 
1580  CALL profiler_set_stop_time('equilibrium_get_gp_ti_pp', &
1581  & start_time)
1582 
1583  END FUNCTION
1584 
1585 !-------------------------------------------------------------------------------
1598 !-------------------------------------------------------------------------------
1599  FUNCTION equilibrium_get_ti_cart(this, x_cart)
1601  IMPLICIT NONE
1602 
1603 ! Declare Arguments
1604  REAL (rprec) :: equilibrium_get_ti_cart
1605  TYPE (equilibrium_class), INTENT(in) :: this
1606  REAL (rprec), DIMENSION(3), INTENT(in) :: x_cart
1607 
1608 ! local variables
1609  REAL (rprec) :: start_time
1610 
1611 ! Start of executable code
1612  start_time = profiler_get_start_time()
1613 
1614  SELECT CASE (this%type)
1615 
1616  CASE (equilibrium_vmec_type)
1617  equilibrium_get_ti_cart = vmec_get_ti(this%vmec, x_cart)
1618 
1620  equilibrium_get_ti_cart = siesta_get_ti(this%siesta, x_cart)
1621 
1622  CASE DEFAULT
1624 
1625  END SELECT
1626 
1627  CALL profiler_set_stop_time('equilibrium_get_ti_cart', &
1628  & start_time)
1629 
1630  END FUNCTION
1631 
1632 !-------------------------------------------------------------------------------
1645 !-------------------------------------------------------------------------------
1646  FUNCTION equilibrium_get_ti_radial(this, r)
1648  IMPLICIT NONE
1649 
1650 ! Declare Arguments
1651  REAL (rprec) :: equilibrium_get_ti_radial
1652  TYPE (equilibrium_class), INTENT(in) :: this
1653  REAL (rprec), INTENT(in) :: r
1654 
1655 ! local variables
1656  REAL (rprec) :: start_time
1657 
1658 ! Start of executable code
1659  start_time = profiler_get_start_time()
1660 
1661  SELECT CASE (this%type)
1662 
1663  CASE (equilibrium_vmec_type)
1664  equilibrium_get_ti_radial = vmec_get_ti(this%vmec, r)
1665 
1667  equilibrium_get_ti_radial = siesta_get_ti(this%siesta, r)
1668 
1669  CASE DEFAULT
1671 
1672  END SELECT
1673 
1674  CALL profiler_set_stop_time('equilibrium_get_ti_radial', &
1675  & start_time)
1676 
1677  END FUNCTION
1678 
1679 !-------------------------------------------------------------------------------
1691 !-------------------------------------------------------------------------------
1692  FUNCTION equilibrium_get_ze_cart(this, x_cart)
1694  IMPLICIT NONE
1695 
1696 ! Declare Arguments
1697  REAL (rprec) :: equilibrium_get_ze_cart
1698  TYPE (equilibrium_class), INTENT(in) :: this
1699  REAL (rprec), DIMENSION(3), INTENT(in) :: x_cart
1700 
1701 ! local variables
1702  REAL (rprec) :: start_time
1703 
1704 ! Start of executable code
1705  start_time = profiler_get_start_time()
1706 
1707  SELECT CASE (this%type)
1708 
1709  CASE (equilibrium_vmec_type)
1710  equilibrium_get_ze_cart = vmec_get_ze(this%vmec, x_cart)
1711 
1712 ! @TODO add ze to siesta
1713 ! CASE (equilibrium_siesta_type)
1714 ! equilibrium_get_ze_cart = siesta_get_ze(this%siesta, x_cart)
1715 
1716  CASE DEFAULT
1718 
1719  END SELECT
1720 
1721  CALL profiler_set_stop_time('equilibrium_get_ze_cart', &
1722  & start_time)
1723 
1724  END FUNCTION
1725 
1726 !-------------------------------------------------------------------------------
1738 !-------------------------------------------------------------------------------
1739  FUNCTION equilibrium_get_ze_radial(this, r)
1741  IMPLICIT NONE
1742 
1743 ! Declare Arguments
1744  REAL (rprec) :: equilibrium_get_ze_radial
1745  TYPE (equilibrium_class), INTENT(in) :: this
1746  REAL (rprec), INTENT(in) :: r
1747 
1748 ! local variables
1749  REAL (rprec) :: start_time
1750 
1751 ! Start of executable code
1752  start_time = profiler_get_start_time()
1753 
1754  SELECT CASE (this%type)
1755 
1756  CASE (equilibrium_vmec_type)
1757  equilibrium_get_ze_radial = vmec_get_ze(this%vmec, r)
1758 
1759 ! @TODO add ze to siesta
1760 ! CASE (equilibrium_siesta_type)
1761 ! equilibrium_get_ze_radial = siesta_get_ze(this%siesta, r)
1762 
1763  CASE DEFAULT
1765 
1766  END SELECT
1767 
1768  CALL profiler_set_stop_time('equilibrium_get_ze_radial', &
1769  & start_time)
1770 
1771  END FUNCTION
1772 
1773 !-------------------------------------------------------------------------------
1784 !-------------------------------------------------------------------------------
1785  FUNCTION equilibrium_get_gp_sxrem_num_hyper_param(this, index)
1787  IMPLICIT NONE
1788 
1789 ! Declare Arguments
1791  TYPE (equilibrium_class), INTENT(in) :: this
1792  INTEGER, INTENT(in) :: index
1793 
1794 ! local variables
1795  REAL (rprec) :: start_time
1796 
1797 ! Start of executable code
1798  start_time = profiler_get_start_time()
1799 
1800  SELECT CASE (this%type)
1801 
1802  CASE (equilibrium_vmec_type)
1804  & vmec_get_gp_sxrem_num_hyper_param(this%vmec, index)
1805 
1808  & siesta_get_gp_sxrem_num_hyper_param(this%siesta, index)
1809 
1810  CASE DEFAULT
1812 
1813  END SELECT
1814 
1815  CALL profiler_set_stop_time( &
1816  & 'equilibrium_get_gp_sxrem_num_hyper_param', start_time)
1817 
1818  END FUNCTION
1819 
1820 !-------------------------------------------------------------------------------
1831 !-------------------------------------------------------------------------------
1832  FUNCTION equilibrium_get_sxrem_af(this, index)
1834  IMPLICIT NONE
1835 
1836 ! Declare Arguments
1837  REAL (rprec), DIMENSION(:), POINTER :: equilibrium_get_sxrem_af
1838  TYPE (equilibrium_class), INTENT(in) :: this
1839  INTEGER, INTENT(in) :: index
1840 
1841 ! local variables
1842  REAL (rprec) :: start_time
1843 
1844 ! Start of executable code
1845  start_time = profiler_get_start_time()
1846 
1847  SELECT CASE (this%type)
1848 
1849  CASE (equilibrium_vmec_type)
1850  equilibrium_get_sxrem_af => vmec_get_sxrem_af(this%vmec, &
1851  & index)
1852 
1855  & index)
1856 
1857  CASE DEFAULT
1858  equilibrium_get_sxrem_af => null()
1859 
1860  END SELECT
1861 
1862  END FUNCTION
1863 
1864 !-------------------------------------------------------------------------------
1879 !-------------------------------------------------------------------------------
1880  FUNCTION equilibrium_get_gp_sxrem_ij(this, i, j, index)
1882  IMPLICIT NONE
1883 
1884 ! Declare Arguments
1885  REAL (rprec) :: equilibrium_get_gp_sxrem_ij
1886  TYPE (equilibrium_class), INTENT(in) :: this
1887  INTEGER, INTENT(in) :: i
1888  INTEGER, INTENT(in) :: j
1889  INTEGER, INTENT(in) :: index
1890 
1891 ! local variables
1892  REAL (rprec) :: start_time
1893 
1894 ! Start of executable code
1895  start_time = profiler_get_start_time()
1896 
1897  SELECT CASE (this%type)
1898 
1899  CASE (equilibrium_vmec_type)
1901  & vmec_get_gp_sxrem(this%vmec, i, j, index)
1902 
1905  & siesta_get_gp_sxrem(this%siesta, i, j, index)
1906 
1907  CASE DEFAULT
1909 
1910  END SELECT
1911 
1912  CALL profiler_set_stop_time('equilibrium_get_gp_sxrem_ij', &
1913  & start_time)
1914 
1915  END FUNCTION
1916 
1917 !-------------------------------------------------------------------------------
1933 !-------------------------------------------------------------------------------
1934  FUNCTION equilibrium_get_gp_sxrem_pi(this, x_cart, i, index)
1936  IMPLICIT NONE
1937 
1938 ! Declare Arguments
1939  REAL (rprec) :: equilibrium_get_gp_sxrem_pi
1940  TYPE (equilibrium_class), INTENT(in) :: this
1941  REAL (rprec), DIMENSION(3), INTENT(in) :: x_cart
1942  INTEGER, INTENT(in) :: i
1943  INTEGER, INTENT(in) :: index
1944 
1945 ! local variables
1946  REAL (rprec) :: start_time
1947 
1948 ! Start of executable code
1949  start_time = profiler_get_start_time()
1950 
1951  SELECT CASE (this%type)
1952 
1953  CASE (equilibrium_vmec_type)
1955  & vmec_get_gp_sxrem(this%vmec, x_cart, i, index)
1956 
1959  & siesta_get_gp_sxrem(this%siesta, x_cart, i, index)
1960 
1961  CASE DEFAULT
1963 
1964  END SELECT
1965 
1966  CALL profiler_set_stop_time('equilibrium_get_gp_sxrem_pi', &
1967  & start_time)
1968 
1969  END FUNCTION
1970 
1971 !-------------------------------------------------------------------------------
1987 !-------------------------------------------------------------------------------
1988  FUNCTION equilibrium_get_gp_sxrem_pp(this, x_cart, y_cart, index)
1990  IMPLICIT NONE
1991 
1992 ! Declare Arguments
1993  REAL (rprec) :: equilibrium_get_gp_sxrem_pp
1994  TYPE (equilibrium_class), INTENT(in) :: this
1995  REAL (rprec), DIMENSION(3), INTENT(in) :: x_cart
1996  REAL (rprec), DIMENSION(3), INTENT(in) :: y_cart
1997  INTEGER, INTENT(in) :: index
1998 
1999 ! local variables
2000  REAL (rprec) :: start_time
2001 
2002 ! Start of executable code
2003  start_time = profiler_get_start_time()
2004 
2005  SELECT CASE (this%type)
2006 
2007  CASE (equilibrium_vmec_type)
2009  & vmec_get_gp_sxrem(this%vmec, x_cart, y_cart, index)
2010 
2013  & siesta_get_gp_sxrem(this%siesta, x_cart, y_cart, index)
2014 
2015  CASE DEFAULT
2017 
2018  END SELECT
2019 
2020  CALL profiler_set_stop_time('equilibrium_get_gp_sxrem_pp', &
2021  & start_time)
2022 
2023  END FUNCTION
2024 
2025 !-------------------------------------------------------------------------------
2039 !-------------------------------------------------------------------------------
2040  FUNCTION equilibrium_get_sxrem_cart(this, x_cart, index)
2042  IMPLICIT NONE
2043 
2044 ! Declare Arguments
2045  REAL (rprec) :: equilibrium_get_sxrem_cart
2046  TYPE (equilibrium_class), INTENT(in) :: this
2047  REAL (rprec), DIMENSION(3), INTENT(in) :: x_cart
2048  INTEGER, INTENT(in) :: index
2049 
2050 ! local variables
2051  REAL (rprec) :: start_time
2052 
2053 ! Start of executable code
2054  start_time = profiler_get_start_time()
2055 
2056  SELECT CASE (this%type)
2057 
2058  CASE (equilibrium_vmec_type)
2059  equilibrium_get_sxrem_cart = vmec_get_sxrem(this%vmec, &
2060  & x_cart, index)
2061 
2064  & x_cart, index)
2065 
2066  CASE DEFAULT
2068 
2069  END SELECT
2070 
2071  CALL profiler_set_stop_time('equilibrium_get_sxrem_cart', &
2072  & start_time)
2073 
2074  END FUNCTION
2075 
2076 !-------------------------------------------------------------------------------
2090 !-------------------------------------------------------------------------------
2091  FUNCTION equilibrium_get_sxrem_radial(this, r, index)
2093  IMPLICIT NONE
2094 
2095 ! Declare Arguments
2096  REAL (rprec) :: equilibrium_get_sxrem_radial
2097  TYPE (equilibrium_class), INTENT(in) :: this
2098  REAL (rprec), INTENT(in) :: r
2099  INTEGER, INTENT(in) :: index
2100 
2101 ! local variables
2102  REAL (rprec) :: start_time
2103 
2104 ! Start of executable code
2105  start_time = profiler_get_start_time()
2106 
2107  SELECT CASE (this%type)
2108 
2109  CASE (equilibrium_vmec_type)
2110  equilibrium_get_sxrem_radial = vmec_get_sxrem(this%vmec, &
2111  & r, index)
2112 
2115  & r, index)
2116 
2117  CASE DEFAULT
2119 
2120  END SELECT
2121 
2122  CALL profiler_set_stop_time('equilibrium_get_sxrem_radial', &
2123  & start_time)
2124 
2125  END FUNCTION
2126 
2127 !-------------------------------------------------------------------------------
2140 !-------------------------------------------------------------------------------
2141  FUNCTION equilibrium_get_p_cart(this, x_cart)
2143  IMPLICIT NONE
2144 
2145 ! Declare Arguments
2146  REAL (rprec) :: equilibrium_get_p_cart
2147  TYPE (equilibrium_class), INTENT(in) :: this
2148  REAL (rprec), DIMENSION(3), INTENT(in) :: x_cart
2149 
2150 ! local variables
2151  REAL (rprec) :: start_time
2152 
2153 ! Start of executable code
2154  start_time = profiler_get_start_time()
2155 
2156  SELECT CASE (this%type)
2157 
2158  CASE (equilibrium_vmec_type)
2159  equilibrium_get_p_cart = vmec_get_p(this%vmec, x_cart)
2160 
2162  equilibrium_get_p_cart = siesta_get_p(this%siesta, x_cart, &
2163  & .false.)
2164 
2165  CASE DEFAULT
2167 
2168  END SELECT
2169 
2170  CALL profiler_set_stop_time('equilibrium_get_p_cart', start_time)
2171 
2172  END FUNCTION
2173 
2174 !-------------------------------------------------------------------------------
2187 !-------------------------------------------------------------------------------
2188  FUNCTION equilibrium_get_p_radial(this, r)
2190  IMPLICIT NONE
2191 
2192 ! Declare Arguments
2193  REAL (rprec) :: equilibrium_get_p_radial
2194  TYPE (equilibrium_class), INTENT(in) :: this
2195  REAL (rprec), INTENT(in) :: r
2196 
2197 ! local variables
2198  REAL (rprec) :: start_time
2199 
2200 ! Start of executable code
2201  start_time = profiler_get_start_time()
2202 
2203  SELECT CASE (this%type)
2204 
2205  CASE (equilibrium_vmec_type)
2206  equilibrium_get_p_radial = vmec_get_p(this%vmec, r)
2207 
2209  equilibrium_get_p_radial = siesta_get_p(this%siesta, r, &
2210  & .false.)
2211 
2212  CASE DEFAULT
2214 
2215  END SELECT
2216 
2217  CALL profiler_set_stop_time('equilibrium_get_p_radial', &
2218  & start_time)
2219 
2220  END FUNCTION
2221 
2222 !-------------------------------------------------------------------------------
2236 !-------------------------------------------------------------------------------
2237  FUNCTION equilibrium_get_b_vec(this, x_cart, cyl)
2239  IMPLICIT NONE
2240 
2241 ! Declare Arguments
2242  REAL (rprec), DIMENSION(3) :: equilibrium_get_b_vec
2243  TYPE (equilibrium_class), INTENT(in) :: this
2244  REAL (rprec), DIMENSION(3), INTENT(in) :: x_cart
2245  LOGICAL, INTENT(in) :: cyl
2246 
2247 ! local variables
2248  REAL (rprec) :: start_time
2249 
2250 ! Start of executable code
2251  start_time = profiler_get_start_time()
2252 
2253  SELECT CASE (this%type)
2254 
2255  CASE (equilibrium_vmec_type)
2256  equilibrium_get_b_vec = vmec_get_b_vec(this%vmec, x_cart, &
2257  & cyl)
2258 
2260  equilibrium_get_b_vec = vacuum_get_b_vec(this%vacuum, &
2261  & x_cart, cyl)
2262 
2264  equilibrium_get_b_vec = siesta_get_b_vec(this%siesta, &
2265  & x_cart, cyl)
2266 
2267  CASE DEFAULT
2268  equilibrium_get_b_vec = 0.0
2269 
2270  END SELECT
2271 
2272  CALL profiler_set_stop_time('equilibrium_get_B_vec', start_time)
2273 
2274  END FUNCTION
2275 
2276 !-------------------------------------------------------------------------------
2290 !-------------------------------------------------------------------------------
2291  FUNCTION equilibrium_get_int_b_dphi(this, r, theta)
2293  IMPLICIT NONE
2294 
2295 ! Declare Arguments
2296  REAL (rprec) :: equilibrium_get_int_b_dphi
2297  TYPE (equilibrium_class), INTENT(in) :: this
2298  REAL (rprec), INTENT(in) :: r
2299  REAL (rprec), INTENT(in) :: theta
2300 
2301 ! local variables
2302  REAL (rprec) :: start_time
2303 
2304 ! Start of executable code
2305  start_time = profiler_get_start_time()
2306 
2307  SELECT CASE (this%type)
2308 
2309  CASE (equilibrium_vmec_type)
2311  & vmec_get_int_b_dphi(this%vmec, r, theta)
2312 
2315  & vacuum_get_int_b_dphi(this%vacuum, r, theta)
2316 
2319  & siesta_get_int_b_dphi(this%siesta, r, theta)
2320 
2321  CASE DEFAULT
2323 
2324  END SELECT
2325 
2326  CALL profiler_set_stop_time('equilibrium_get_Int_B_dphi', &
2327  & start_time)
2328 
2329  END FUNCTION
2330 
2331 !-------------------------------------------------------------------------------
2352 !-------------------------------------------------------------------------------
2353  FUNCTION equilibrium_get_plasma_edge(this, phi, r, z)
2355  IMPLICIT NONE
2356 
2357 ! Declare Arguments
2358  INTEGER :: equilibrium_get_plasma_edge
2359  TYPE (equilibrium_class), INTENT(in) :: this
2360  REAL (rprec), INTENT (in) :: phi
2361  REAL (rprec), DIMENSION(:), POINTER :: r
2362  REAL (rprec), DIMENSION(:), POINTER :: z
2363 
2364 ! local variables
2365  REAL (rprec) :: start_time
2366 
2367 ! Start of executable code
2368  start_time = profiler_get_start_time()
2369 
2370  SELECT CASE (this%type)
2371 
2372  CASE (equilibrium_vmec_type)
2374  & vmec_get_plasma_edge(this%vmec, phi, r, z)
2375 
2378  & siesta_get_plasma_edge(this%siesta, phi, r, z)
2379 
2380  CASE DEFAULT
2382  r => null()
2383  z => null()
2384 
2385  END SELECT
2386 
2387  CALL profiler_set_stop_time('equilibrium_get_plasma_edge', &
2388  & start_time)
2389 
2390  END FUNCTION
2391 
2392 !-------------------------------------------------------------------------------
2402 !-------------------------------------------------------------------------------
2405  IMPLICIT NONE
2406 
2407 ! Declare Arguments
2408  REAL (rprec), DIMENSION(:,:,:), POINTER :: &
2410  TYPE (equilibrium_class), INTENT(in) :: this
2411 
2412 ! local variables
2413  REAL (rprec) :: start_time
2414 
2415 ! Start of executable code
2416  start_time = profiler_get_start_time()
2417 
2418  SELECT CASE (this%type)
2419 
2420  CASE (equilibrium_vmec_type)
2422  & vmec_get_magnetic_volume_rgrid(this%vmec)
2423 
2426  & siesta_get_magnetic_volume_rgrid(this%siesta)
2427 
2428  CASE DEFAULT
2430 
2431  END SELECT
2432 
2433  CALL profiler_set_stop_time( &
2434  & 'equilibrium_get_magnetic_volume_rgrid', start_time)
2435 
2436  END FUNCTION
2437 
2438 !-------------------------------------------------------------------------------
2448 !-------------------------------------------------------------------------------
2451  IMPLICIT NONE
2452 
2453 ! Declare Arguments
2454  REAL (rprec), DIMENSION(:,:,:), POINTER :: &
2456  TYPE (equilibrium_class), INTENT(in) :: this
2457 
2458 ! local variables
2459  REAL (rprec) :: start_time
2460 
2461 ! Start of executable code
2462  start_time = profiler_get_start_time()
2463 
2464  SELECT CASE (this%type)
2465 
2466  CASE (equilibrium_vmec_type)
2468  & vmec_get_magnetic_volume_zgrid(this%vmec)
2469 
2472  & siesta_get_magnetic_volume_zgrid(this%siesta)
2473 
2474  CASE DEFAULT
2476 
2477  END SELECT
2478 
2479  CALL profiler_set_stop_time( &
2480  & 'equilibrium_get_magnetic_volume_zgrid', start_time)
2481 
2482  END FUNCTION
2483 
2484 !-------------------------------------------------------------------------------
2494 !-------------------------------------------------------------------------------
2497  IMPLICIT NONE
2498 
2499 ! Declare Arguments
2500  REAL (rprec), DIMENSION(:,:,:), POINTER :: &
2502  TYPE (equilibrium_class), INTENT(in) :: this
2503 
2504 ! local variables
2505  REAL (rprec) :: start_time
2506 
2507 ! Start of executable code
2508  start_time = profiler_get_start_time()
2509 
2510  SELECT CASE (this%type)
2511 
2512  CASE (equilibrium_vmec_type)
2514  & vmec_get_magnetic_volume_jrgrid(this%vmec)
2515 
2518  & siesta_get_magnetic_volume_jrgrid(this%siesta)
2519 
2520  CASE DEFAULT
2522 
2523  END SELECT
2524 
2525  CALL profiler_set_stop_time( &
2526  & 'equilibrium_get_magnetic_volume_jrgrid', start_time)
2527 
2528  END FUNCTION
2529 
2530 !-------------------------------------------------------------------------------
2540 !-------------------------------------------------------------------------------
2543  IMPLICIT NONE
2544 
2545 ! Declare Arguments
2546  REAL (rprec), DIMENSION(:,:,:), POINTER :: &
2548  TYPE (equilibrium_class), INTENT(in) :: this
2549 
2550 ! local variables
2551  REAL (rprec) :: start_time
2552 
2553 ! Start of executable code
2554  start_time = profiler_get_start_time()
2555 
2556  SELECT CASE (this%type)
2557 
2558  CASE (equilibrium_vmec_type)
2560  & vmec_get_magnetic_volume_jphigrid(this%vmec)
2561 
2565 
2566  CASE DEFAULT
2568 
2569  END SELECT
2570 
2571  CALL profiler_set_stop_time( &
2572  & 'equilibrium_get_magnetic_volume_jphigrid', start_time)
2573 
2574  END FUNCTION
2575 
2576 !-------------------------------------------------------------------------------
2586 !-------------------------------------------------------------------------------
2589  IMPLICIT NONE
2590 
2591 ! Declare Arguments
2592  REAL (rprec), DIMENSION(:,:,:), POINTER :: &
2594  TYPE (equilibrium_class), INTENT(in) :: this
2595 
2596 ! local variables
2597  REAL (rprec) :: start_time
2598 
2599 ! Start of executable code
2600  start_time = profiler_get_start_time()
2601 
2602  SELECT CASE (this%type)
2603 
2604  CASE (equilibrium_vmec_type)
2606  & vmec_get_magnetic_volume_jzgrid(this%vmec)
2607 
2610  & siesta_get_magnetic_volume_jzgrid(this%siesta)
2611 
2612  CASE DEFAULT
2614 
2615  END SELECT
2616 
2617  CALL profiler_set_stop_time( &
2618  & 'equilibrium_get_magnetic_volume_jzgrid', start_time)
2619 
2620  END FUNCTION
2621 
2622 !-------------------------------------------------------------------------------
2632 !-------------------------------------------------------------------------------
2633  FUNCTION equilibrium_get_volume_int_element(this)
2635  IMPLICIT NONE
2636 
2637 ! Declare Arguments
2638  REAL (rprec) :: equilibrium_get_volume_int_element
2639  TYPE (equilibrium_class), INTENT(in) :: this
2640 
2641 ! local variables
2642  REAL (rprec) :: start_time
2643 
2644 ! Start of executable code
2645  start_time = profiler_get_start_time()
2646 
2647  SELECT CASE (this%type)
2648 
2649  CASE (equilibrium_vmec_type)
2651  & vmec_get_volume_int_element(this%vmec)
2652 
2655  & siesta_get_volume_int_element(this%siesta)
2656 
2657  CASE DEFAULT
2659 
2660  END SELECT
2661 
2662  CALL profiler_set_stop_time('equilibrium_get_volume_int_element', &
2663  & start_time)
2664 
2665  END FUNCTION
2666 
2667 !-------------------------------------------------------------------------------
2677 !-------------------------------------------------------------------------------
2678  FUNCTION equilibrium_get_con_surface_krgrid(this)
2680  IMPLICIT NONE
2681 
2682 ! Declare Arguments
2683  REAL (rprec), DIMENSION(:,:), POINTER :: &
2685  TYPE (equilibrium_class), INTENT(in) :: this
2686 
2687 ! local variables
2688  REAL (rprec) :: start_time
2689 
2690 ! Start of executable code
2691  start_time = profiler_get_start_time()
2692 
2693  SELECT CASE (this%type)
2694 
2695  CASE (equilibrium_vmec_type)
2697  & vmec_get_con_surface_krgrid(this%vmec)
2698 
2701  & siesta_get_con_surface_krgrid(this%siesta)
2702 
2703  CASE DEFAULT
2705 
2706  END SELECT
2707 
2708  CALL profiler_set_stop_time( &
2709  & 'equilibrium_get_con_surface_krgrid', start_time)
2710 
2711  END FUNCTION
2712 
2713 !-------------------------------------------------------------------------------
2723 !-------------------------------------------------------------------------------
2726  IMPLICIT NONE
2727 
2728 ! Declare Arguments
2729  REAL (rprec), DIMENSION(:,:), POINTER :: &
2731  TYPE (equilibrium_class), INTENT(in) :: this
2732 
2733 ! local variables
2734  REAL (rprec) :: start_time
2735 
2736 ! Start of executable code
2737  start_time = profiler_get_start_time()
2738 
2739  SELECT CASE (this%type)
2740 
2741  CASE (equilibrium_vmec_type)
2743  & vmec_get_con_surface_kphigrid(this%vmec)
2744 
2747  & siesta_get_con_surface_kphigrid(this%siesta)
2748 
2749  CASE DEFAULT
2751 
2752  END SELECT
2753 
2754  CALL profiler_set_stop_time( &
2755  & 'equilibrium_get_con_surface_kphigrid', start_time)
2756 
2757  END FUNCTION
2758 
2759 !-------------------------------------------------------------------------------
2769 !-------------------------------------------------------------------------------
2770  FUNCTION equilibrium_get_con_surface_kzgrid(this)
2772  IMPLICIT NONE
2773 
2774 ! Declare Arguments
2775  REAL (rprec), DIMENSION(:,:), POINTER :: &
2777  TYPE (equilibrium_class), INTENT(in) :: this
2778 
2779 ! local variables
2780  REAL (rprec) :: start_time
2781 
2782 ! Start of executable code
2783  start_time = profiler_get_start_time()
2784 
2785  SELECT CASE (this%type)
2786 
2787  CASE (equilibrium_vmec_type)
2789  & vmec_get_con_surface_kzgrid(this%vmec)
2790 
2793  & siesta_get_con_surface_kzgrid(this%siesta)
2794 
2795  CASE DEFAULT
2797 
2798  END SELECT
2799 
2800  CALL profiler_set_stop_time( &
2801  & 'equilibrium_get_con_surface_kzgrid', start_time)
2802 
2803  END FUNCTION
2804 
2805 !-------------------------------------------------------------------------------
2815 !-------------------------------------------------------------------------------
2816  FUNCTION equilibrium_get_area_int_element(this)
2818  IMPLICIT NONE
2819 
2820 ! Declare Arguments
2821  REAL (rprec) :: equilibrium_get_area_int_element
2822  TYPE (equilibrium_class), INTENT(in) :: this
2823 
2824 ! local variables
2825  REAL (rprec) :: start_time
2826 
2827 ! Start of executable code
2828  start_time = profiler_get_start_time()
2829 
2830  SELECT CASE (this%type)
2831 
2832  CASE (equilibrium_vmec_type)
2834  & vmec_get_area_int_element(this%vmec)
2835 
2838  & siesta_get_area_int_element(this%siesta)
2839 
2840  CASE DEFAULT
2842 
2843  END SELECT
2844 
2845  CALL profiler_set_stop_time('equilibrium_get_area_int_element', &
2846  & start_time)
2847 
2848  END FUNCTION
2849 
2850 !-------------------------------------------------------------------------------
2870 !-------------------------------------------------------------------------------
2871  FUNCTION equilibrium_get_ext_currents(this, num_currents, &
2872  & scale_currents)
2874  IMPLICIT NONE
2875 
2876 ! Declare Arguments
2877  REAL (rprec), DIMENSION(:), POINTER :: &
2879  TYPE (equilibrium_class), INTENT(in) :: this
2880  INTEGER :: num_currents
2881  LOGICAL, INTENT(out) :: scale_currents
2882 
2883 ! local variables
2884  REAL (rprec) :: start_time
2885 
2886 ! Start of executable code
2887  start_time = profiler_get_start_time()
2888 
2889  SELECT CASE (this%type)
2890 
2891  CASE (equilibrium_vmec_type)
2893  & vmec_get_ext_currents(this%vmec, num_currents, &
2894  & scale_currents)
2895 
2897 ! Vacuum equilibira always return currents. No need to force it.
2899  & vacuum_get_ext_currents(this%vacuum, scale_currents)
2900 
2903  & siesta_get_ext_currents(this%siesta, num_currents, &
2904  & scale_currents)
2905 
2906  CASE DEFAULT
2907  scale_currents = .false.
2909 
2910  END SELECT
2911 
2912  CALL profiler_set_stop_time('equilibrium_get_ext_currents', &
2913  & start_time)
2914 
2915  END FUNCTION
2916 
2917 !-------------------------------------------------------------------------------
2930 !-------------------------------------------------------------------------------
2931  FUNCTION equilibrium_get_ext_b_plasma(this, position, axi_only)
2933  IMPLICIT NONE
2934 
2935 ! Declare Arguments
2936  REAL (rprec), DIMENSION(3) :: equilibrium_get_ext_b_plasma
2937  TYPE (equilibrium_class), INTENT(in) :: this
2938  REAL (rprec), DIMENSION(3) :: position
2939  LOGICAL, INTENT(in) :: axi_only
2940 
2941 ! local variables
2942  REAL (rprec) :: start_time
2943 
2944 ! Start of executable code
2945  start_time = profiler_get_start_time()
2946 
2947  SELECT CASE (this%type)
2948 
2949  CASE (equilibrium_vmec_type)
2951  & vmec_get_ext_b_plasma(this%vmec, position, axi_only)
2952 
2955  & siesta_get_ext_b_plasma(this%siesta, position, axi_only)
2956 
2957  CASE DEFAULT
2959 
2960  END SELECT
2961 
2962  CALL profiler_set_stop_time('equilibrium_get_ext_b_plasma', &
2963  & start_time)
2964 
2965  END FUNCTION
2966 
2967 !-------------------------------------------------------------------------------
2977 !-------------------------------------------------------------------------------
2978  FUNCTION equilibrium_get_grid_size(this)
2980  IMPLICIT NONE
2981 
2982 ! Declare Arguments
2983  INTEGER :: equilibrium_get_grid_size
2984  TYPE (equilibrium_class), INTENT(in) :: this
2985 
2986 ! local variables
2987  REAL (rprec) :: start_time
2988 
2989 ! Start of executable code
2990  start_time = profiler_get_start_time()
2991 
2992  SELECT CASE (this%type)
2993 
2994  CASE (equilibrium_vmec_type)
2995  equilibrium_get_grid_size = vmec_get_grid_size(this%vmec)
2996 
2999  & siesta_get_grid_size(this%siesta)
3000 
3001  CASE DEFAULT
3003 
3004  END SELECT
3005 
3006  CALL profiler_set_stop_time('equilibrium_get_grid_size', &
3007  & start_time)
3008 
3009  END FUNCTION
3010 
3011 !-------------------------------------------------------------------------------
3021 !-------------------------------------------------------------------------------
3022  FUNCTION equilibrium_get_grid_start(this)
3024  IMPLICIT NONE
3025 
3026 ! Declare Arguments
3027  REAl (rprec) :: equilibrium_get_grid_start
3028  TYPE (equilibrium_class), INTENT(in) :: this
3029 
3030 ! local variables
3031  REAL (rprec) :: start_time
3032 
3033 ! Start of executable code
3034  start_time = profiler_get_start_time()
3035 
3036  SELECT CASE (this%type)
3037 
3038  CASE (equilibrium_vmec_type)
3039  equilibrium_get_grid_start = vmec_get_grid_start(this%vmec)
3040 
3043  & siesta_get_grid_start(this%siesta)
3044 
3045  CASE DEFAULT
3047 
3048  END SELECT
3049 
3050  CALL profiler_set_stop_time('equilibrium_get_grid_start', &
3051  & start_time)
3052 
3053  END FUNCTION
3054 
3055 !-------------------------------------------------------------------------------
3065 !-------------------------------------------------------------------------------
3066  FUNCTION equilibrium_get_grid_stop(this)
3068  IMPLICIT NONE
3069 
3070 ! Declare Arguments
3071  REAl (rprec) :: equilibrium_get_grid_stop
3072  TYPE (equilibrium_class), INTENT(in) :: this
3073 
3074 ! local variables
3075  REAL (rprec) :: start_time
3076 
3077 ! Start of executable code
3078  start_time = profiler_get_start_time()
3079 
3080  SELECT CASE (this%type)
3081 
3082  CASE (equilibrium_vmec_type)
3083  equilibrium_get_grid_stop = vmec_get_grid_stop(this%vmec)
3084 
3087  & siesta_get_grid_stop(this%siesta)
3088 
3089  CASE DEFAULT
3091 
3092  END SELECT
3093 
3094  CALL profiler_set_stop_time('equilibrium_get_grid_stop', &
3095  & start_time)
3096 
3097  END FUNCTION
3098 
3099 !*******************************************************************************
3100 ! QUERY SUBROUTINES
3101 !*******************************************************************************
3102 !-------------------------------------------------------------------------------
3113 !-------------------------------------------------------------------------------
3114  FUNCTION equilibrium_is_scaler_value(this, id)
3116  IMPLICIT NONE
3117 
3118 ! Declare Arguments
3119  LOGICAL :: equilibrium_is_scaler_value
3120  TYPE (equilibrium_class), INTENT(in) :: this
3121  INTEGER, INTENT(in) :: id
3122 
3123 ! local variables
3124  REAL (rprec) :: start_time
3125 
3126 ! Start of executable code
3127  start_time = profiler_get_start_time()
3128 
3129  SELECT CASE (this%type)
3130 
3131  CASE (equilibrium_vmec_type)
3133  & vmec_is_scaler_value(this%vmec, id)
3134 
3137  & siesta_is_scaler_value(this%siesta, id)
3138 
3139  CASE DEFAULT
3140  equilibrium_is_scaler_value = .false.
3141 
3142  END SELECT
3143 
3144  CALL profiler_set_stop_time('equilibrium_is_scaler_value', &
3145  & start_time)
3146 
3147  END FUNCTION
3148 
3149 !-------------------------------------------------------------------------------
3161 !-------------------------------------------------------------------------------
3162  FUNCTION equilibrium_is_1d_array(this, id)
3164  IMPLICIT NONE
3165 
3166 ! Declare Arguments
3167  LOGICAL :: equilibrium_is_1d_array
3168  TYPE (equilibrium_class), INTENT(in) :: this
3169  INTEGER, INTENT(in) :: id
3170 
3171 ! local variables
3172  REAL (rprec) :: start_time
3173 
3174 ! Start of executable code
3175  start_time = profiler_get_start_time()
3176 
3177  SELECT CASE (this%type)
3178 
3179  CASE (equilibrium_vmec_type)
3180  equilibrium_is_1d_array = vmec_is_1d_array(this%vmec, id)
3181 
3184  & id)
3185 
3188  & id)
3189 
3190  CASE DEFAULT
3191  equilibrium_is_1d_array = .false.
3192 
3193  END SELECT
3194 
3195  CALL profiler_set_stop_time('equilibrium_is_1d_array', start_time)
3196 
3197  END FUNCTION
3198 
3199 !-------------------------------------------------------------------------------
3210 !-------------------------------------------------------------------------------
3211  FUNCTION equilibrium_is_2d_array(this, id)
3213  IMPLICIT NONE
3214 
3215 ! Declare Arguments
3216  LOGICAL :: equilibrium_is_2d_array
3217  TYPE (equilibrium_class), INTENT(in) :: this
3218  INTEGER, INTENT(in) :: id
3219 
3220 ! local variables
3221  REAL (rprec) :: start_time
3222 
3223 ! Start of executable code
3224  start_time = profiler_get_start_time()
3225 
3226  SELECT CASE (this%type)
3227 
3228  CASE (equilibrium_vmec_type)
3229  equilibrium_is_2d_array = vmec_is_2d_array(this%vmec, id)
3230 
3233  & id)
3234 
3235  CASE DEFAULT
3236  equilibrium_is_2d_array = .false.
3237 
3238  END SELECT
3239 
3240  CALL profiler_set_stop_time('equilibrium_is_2d_array', start_time)
3241 
3242  END FUNCTION
3243 
3244 !-------------------------------------------------------------------------------
3257 !-------------------------------------------------------------------------------
3258  FUNCTION equilibrium_is_recon_param(this, id)
3260  IMPLICIT NONE
3261 
3262 ! Declare Arguments
3263  LOGICAL :: equilibrium_is_recon_param
3264  TYPE (equilibrium_class), INTENT(in) :: this
3265  INTEGER, INTENT(in) :: id
3266 
3267 ! local variables
3268  REAL (rprec) :: start_time
3269 
3270 ! Start of executable code
3271  start_time = profiler_get_start_time()
3272 
3273  SELECT CASE (this%type)
3274 
3275  CASE (equilibrium_vmec_type)
3277  & vmec_is_recon_param(this%vmec, id)
3278 
3281  & vacuum_is_recon_param(this%vacuum, id)
3282 
3285  & siesta_is_recon_param(this%siesta, id)
3286 
3287  CASE DEFAULT
3288  equilibrium_is_recon_param = .false.
3289 
3290  END SELECT
3291 
3292  CALL profiler_set_stop_time('equilibrium_is_recon_param', &
3293  & start_time)
3294 
3295  END FUNCTION
3296 
3297 !-------------------------------------------------------------------------------
3307 !-------------------------------------------------------------------------------
3308  FUNCTION equilibrium_is_using_point(this)
3310  IMPLICIT NONE
3311 
3312 ! Declare Arguments
3313  LOGICAL :: equilibrium_is_using_point
3314  TYPE (equilibrium_class), INTENT(in) :: this
3315 
3316 ! local variables
3317  REAL (rprec) :: start_time
3318 
3319 ! Start of executable code
3320  start_time = profiler_get_start_time()
3321 
3322  SELECT CASE (this%type)
3323 
3324  CASE (equilibrium_vmec_type)
3325  equilibrium_is_using_point = vmec_is_using_point(this%vmec)
3326 
3329  & siesta_is_using_point(this%siesta)
3330 
3331  CASE DEFAULT
3332  equilibrium_is_using_point = .false.
3333 
3334  END SELECT
3335 
3336  CALL profiler_set_stop_time('equilibrium_is_using_point', &
3337  & start_time)
3338 
3339  END FUNCTION
3340 
3341 !*******************************************************************************
3342 ! UTILITY SUBROUTINES
3343 !*******************************************************************************
3344 !-------------------------------------------------------------------------------
3359 !-------------------------------------------------------------------------------
3360  FUNCTION equilibrium_converge(this, num_iter, iou, eq_comm, &
3361  & state_flags)
3363  IMPLICIT NONE
3364 
3365 ! Declare Arguments
3366  LOGICAL :: equilibrium_converge
3367  TYPE (equilibrium_class), INTENT(inout) :: this
3368  INTEGER, INTENT(inout) :: num_iter
3369  INTEGER, INTENT(in) :: iou
3370  INTEGER, INTENT(in) :: eq_comm
3371  INTEGER, INTENT(in) :: state_flags
3372 
3373 ! local variables
3374  INTEGER :: error
3375  INTEGER :: eq_rank
3376  REAL (rprec) :: start_time
3377 
3378 ! Start of executable code
3379  start_time = profiler_get_start_time()
3380 
3381 #if defined(MPI_OPT)
3382  CALL mpi_bcast(num_iter, 1, mpi_integer, 0, eq_comm, error)
3383 #endif
3384  SELECT CASE (this%type)
3385 
3386  CASE (equilibrium_vmec_type)
3387  equilibrium_converge = vmec_converge(this%vmec, num_iter, &
3388  & iou, eq_comm)
3389 
3391  equilibrium_converge = siesta_converge(this%siesta, &
3392  & num_iter, iou, &
3393  & eq_comm, state_flags)
3394 
3395  CASE DEFAULT
3396  equilibrium_converge = .false.
3397 
3398  END SELECT
3399 
3400  CALL profiler_set_stop_time('equilibrium_converge', start_time)
3401 
3402  END FUNCTION
3403 
3404 !-------------------------------------------------------------------------------
3415 !-------------------------------------------------------------------------------
3416  SUBROUTINE equilibrium_read_vac_file(this, index, eq_comm)
3418  IMPLICIT NONE
3419 
3420 ! Declare Arguments
3421  TYPE (equilibrium_class), INTENT(in) :: this
3422  INTEGER, INTENT(in) :: index
3423  INTEGER, INTENT(in) :: eq_comm
3424 
3425 ! Start of executable code
3426  SELECT CASE (this%type)
3427 
3428  CASE (equilibrium_vmec_type)
3429  CALL vmec_read_vac_file(this%vmec, index, eq_comm)
3430 
3432  CALL siesta_read_vac_file(this%siesta, index, eq_comm)
3433 
3434  END SELECT
3435 
3436  END SUBROUTINE
3437 
3438 !-------------------------------------------------------------------------------
3447 !-------------------------------------------------------------------------------
3448  SUBROUTINE equilibrium_save_state(this)
3450  IMPLICIT NONE
3451 
3452 ! Declare Arguments
3453  TYPE (equilibrium_class), INTENT(inout) :: this
3454 
3455 ! local variables
3456  REAL (rprec) :: start_time
3457 
3458 ! Start of executable code
3459  start_time = profiler_get_start_time()
3460 
3461  SELECT CASE (this%type)
3462 
3463  CASE (equilibrium_vmec_type)
3464  CALL vmec_save_state(this%vmec)
3465 
3467  CALL siesta_save_state(this%siesta)
3468 
3469  END SELECT
3470 
3471  CALL profiler_set_stop_time('equilibrium_save_state', start_time)
3472 
3473  END SUBROUTINE
3474 
3475 !-------------------------------------------------------------------------------
3485 !-------------------------------------------------------------------------------
3486  SUBROUTINE equilibrium_reset_state(this)
3488  IMPLICIT NONE
3489 
3490 ! Declare Arguments
3491  TYPE (equilibrium_class), INTENT(inout) :: this
3492 
3493 ! local variables
3494  REAL (rprec) :: start_time
3495 
3496 ! Start of executable code
3497  start_time = profiler_get_start_time()
3498 
3499  SELECT CASE (this%type)
3500 
3501  CASE (equilibrium_vmec_type)
3502  CALL vmec_reset_state(this%vmec)
3503 
3505  CALL siesta_reset_state(this%siesta)
3506 
3507  END SELECT
3508 
3509  CALL profiler_set_stop_time('equilibrium_reset_state', &
3510  & start_time)
3511 
3512  END SUBROUTINE
3513 
3514 !-------------------------------------------------------------------------------
3525 !-------------------------------------------------------------------------------
3526  SUBROUTINE equilibrium_write(this, iou)
3528  IMPLICIT NONE
3529 
3530 ! Declare Arguments
3531  TYPE (equilibrium_class), INTENT(in) :: this
3532  INTEGER, INTENT(in) :: iou
3533 
3534 ! local variables
3535  REAL (rprec) :: start_time
3536 
3537 ! Start of executable code
3538  start_time = profiler_get_start_time()
3539 
3540  SELECT CASE (this%type)
3541 
3542  CASE (equilibrium_vmec_type)
3543  CALL vmec_write(this%vmec, iou)
3544 
3546  CALL vacuum_write(this%vacuum, iou)
3547 
3549  CALL siesta_write(this%siesta, iou)
3550 
3551  END SELECT
3552 
3553  CALL profiler_set_stop_time('equilibrium_write', start_time)
3554 
3555  END SUBROUTINE
3556 
3557 !-------------------------------------------------------------------------------
3568 !-------------------------------------------------------------------------------
3569  SUBROUTINE equilibrium_write_input(this, current_step)
3571  IMPLICIT NONE
3572 
3573 ! Declare Arguments
3574  TYPE (equilibrium_class), INTENT(in) :: this
3575  INTEGER :: current_step
3576 
3577 ! local variables
3578  REAL (rprec) :: start_time
3579 
3580 ! Start of executable code
3581  start_time = profiler_get_start_time()
3582 
3583  SELECT CASE (this%type)
3584 
3585  CASE (equilibrium_vmec_type)
3586  CALL vmec_write_input(this%vmec, current_step)
3587 
3589  CALL vacuum_write_input(this%vacuum, current_step)
3590 
3592  CALL siesta_write_input(this%siesta, current_step)
3593 
3594  END SELECT
3595 
3596  CALL profiler_set_stop_time('equilibrium_write_input', start_time)
3597 
3598  END SUBROUTINE
3599 
3600 !*******************************************************************************
3601 ! NETCDF SUBROUTINES
3602 !*******************************************************************************
3603 !-------------------------------------------------------------------------------
3615 !-------------------------------------------------------------------------------
3616  SUBROUTINE equilibrium_def_result(this, result_ncid, &
3617  & maxnsetps_dim_id)
3619  IMPLICIT NONE
3620 
3621 ! Declare Arguments
3622  TYPE (equilibrium_class), INTENT(in) :: this
3623  INTEGER, INTENT(in) :: result_ncid
3624  INTEGER, INTENT(in) :: maxnsetps_dim_id
3625 
3626 ! local variables
3627  REAL (rprec) :: start_time
3628 
3629 ! Start of executable code
3630  start_time = profiler_get_start_time()
3631 
3632  SELECT CASE (this%type)
3633 
3634  CASE (equilibrium_vmec_type)
3635  CALL vmec_def_result(this%vmec, result_ncid, &
3636  & maxnsetps_dim_id)
3637 
3639  CALL siesta_def_result(this%siesta, result_ncid, &
3640  & maxnsetps_dim_id)
3641 
3642  END SELECT
3643 
3644  CALL profiler_set_stop_time('equilibrium_def_result', start_time)
3645 
3646  END SUBROUTINE
3647 
3648 !-------------------------------------------------------------------------------
3658 !-------------------------------------------------------------------------------
3659  SUBROUTINE equilibrium_write_init_data(this, result_ncid)
3661  IMPLICIT NONE
3662 
3663 ! Declare Arguments
3664  TYPE (equilibrium_class), INTENT(in) :: this
3665  INTEGER, INTENT(in) :: result_ncid
3666 
3667 ! local variables
3668  REAL (rprec) :: start_time
3669 
3670 ! Start of executable code
3671  start_time = profiler_get_start_time()
3672 
3673  SELECT CASE (this%type)
3674 
3675  CASE (equilibrium_vmec_type)
3676  CALL vmec_write_init_data(this%vmec, result_ncid)
3677 
3679  CALL siesta_write_init_data(this%siesta, result_ncid)
3680 
3681  END SELECT
3682 
3683  CALL profiler_set_stop_time('equilibrium_write_init_data', &
3684  & start_time)
3685 
3686  END SUBROUTINE
3687 
3688 !-------------------------------------------------------------------------------
3699 !-------------------------------------------------------------------------------
3700  SUBROUTINE equilibrium_write_step_data(this, result_ncid, &
3701  & current_step)
3703  IMPLICIT NONE
3704 
3705 ! Declare Arguments
3706  TYPE (equilibrium_class), INTENT(in) :: this
3707  INTEGER, INTENT(in) :: result_ncid
3708  INTEGER, INTENT(in) :: current_step
3709 
3710 ! local variables
3711  REAL (rprec) :: start_time
3712 
3713 ! Start of executable code
3714  start_time = profiler_get_start_time()
3715 
3716  SELECT CASE (this%type)
3717 
3718  CASE (equilibrium_vmec_type)
3719  CALL vmec_write_step_data(this%vmec, result_ncid, &
3720  & current_step)
3721 
3723  CALL siesta_write_step_data(this%siesta, result_ncid, &
3724  & current_step)
3725 
3726  END SELECT
3727 
3728  CALL profiler_set_stop_time('equilibrium_write_step_data', &
3729  & start_time)
3730 
3731  END SUBROUTINE
3732 
3733 !-------------------------------------------------------------------------------
3744 !-------------------------------------------------------------------------------
3745  SUBROUTINE equilibrium_restart(this, result_ncid, current_step)
3747  IMPLICIT NONE
3748 
3749 ! Declare Arguments
3750  TYPE (equilibrium_class), INTENT(in) :: this
3751  INTEGER, INTENT(in) :: result_ncid
3752  INTEGER, INTENT(in) :: current_step
3753 
3754 ! local variables
3755  REAL (rprec) :: start_time
3756 
3757 ! Start of executable code
3758  start_time = profiler_get_start_time()
3759 
3760  SELECT CASE (this%type)
3761 
3762  CASE (equilibrium_vmec_type)
3763  CALL vmec_restart(this%vmec, result_ncid, current_step)
3764 
3766  CALL siesta_restart(this%siesta, result_ncid, current_step)
3767 
3768  END SELECT
3769 
3770  CALL profiler_set_stop_time('equilibrium_restart', start_time)
3771 
3772  END SUBROUTINE
3773 
3774 !*******************************************************************************
3775 ! MPI SUBROUTINES
3776 !*******************************************************************************
3777 !-------------------------------------------------------------------------------
3785 !-------------------------------------------------------------------------------
3786  SUBROUTINE equilibrium_sync_state(this, recon_comm)
3788  IMPLICIT NONE
3789 
3790 ! Declare Arguments
3791  TYPE (equilibrium_class), INTENT(inout) :: this
3792  INTEGER, INTENT(in) :: recon_comm
3793 
3794 #if defined(MPI_OPT)
3795 ! local variables
3796  REAL (rprec) :: start_time
3797  INTEGER :: error
3798 
3799 ! Start of executable code
3800  start_time = profiler_get_start_time()
3801 
3802  SELECT CASE (this%type)
3803 
3804  CASE (equilibrium_vmec_type)
3805  CALL vmec_sync_state(this%vmec, recon_comm)
3806 
3808  CALL siesta_sync_state(this%siesta, recon_comm)
3809 
3810  END SELECT
3811 
3812  CALL profiler_set_stop_time('equilibrium_sync_state', start_time)
3813 
3814 #endif
3815  END SUBROUTINE
3816 
3817 !-------------------------------------------------------------------------------
3826 !-------------------------------------------------------------------------------
3827  SUBROUTINE equilibrium_sync_child(this, index, recon_comm)
3829  IMPLICIT NONE
3830 
3831 ! Declare Arguments
3832  TYPE (equilibrium_class), INTENT(inout) :: this
3833  INTEGER, INTENT(in) :: index
3834  INTEGER, INTENT(in) :: recon_comm
3835 
3836 #if defined(MPI_OPT)
3837 ! local variables
3838  REAL (rprec) :: start_time
3839  INTEGER :: mpi_rank
3840  INTEGER :: error
3841 
3842 ! Start of executable code
3843  start_time = profiler_get_start_time()
3844 
3845  CALL mpi_comm_rank(recon_comm, mpi_rank, error)
3846 
3847  SELECT CASE (this%type)
3848 
3849  CASE (equilibrium_vmec_type)
3850  CALL vmec_sync_child(this%vmec, index, recon_comm)
3851 
3853  CALL siesta_sync_child(this%siesta, index, recon_comm)
3854 
3855  END SELECT
3856 
3857  CALL profiler_set_stop_time('equilibrium_sync_child', start_time)
3858 
3859 #endif
3860  END SUBROUTINE
3861 
3862  END MODULE
equilibrium::equilibrium_get_param_name
character(len=data_name_length) function equilibrium_get_param_name(this, id)
Gets the name of a reconstruction equilibrium parameter.
Definition: equilibrium.f:640
equilibrium::equilibrium_get_area_int_element
real(rprec) function equilibrium_get_area_int_element(this)
Get area integration element.
Definition: equilibrium.f:2817
vacuum_equilibrium::vacuum_is_recon_param
logical function vacuum_is_recon_param(this, id)
Checks if a parameter id is a reconstruction parameter.
Definition: vacuum_equilibrium.f:530
equilibrium::equilibrium_get_gp_ti_num_hyper_param
integer function equilibrium_get_gp_ti_num_hyper_param(this)
Get the number of ion temperature gp kernel hyper parameters.
Definition: equilibrium.f:1357
siesta_equilibrium::siesta_get_con_surface_kphigrid
real(rprec) function, dimension(:,:), pointer siesta_get_con_surface_kphigrid(this)
Get the conducting surface integration kphi grid points.
Definition: siesta_equilibrium.f:3381
equilibrium::equilibrium_def_result
subroutine equilibrium_def_result(this, result_ncid, maxnsetps_dim_id)
Define NetCDF variables for the result file.
Definition: equilibrium.f:3618
siesta_equilibrium::siesta_get_magnetic_volume_rgrid
real(rprec) function, dimension(:,:,:), pointer siesta_get_magnetic_volume_rgrid(this)
Get volume magnetic volume integration radial grid points.
Definition: siesta_equilibrium.f:3163
equilibrium::equilibrium_is_2d_array
logical function equilibrium_is_2d_array(this, id)
Checks if a parameter id is a 2d array.
Definition: equilibrium.f:3212
siesta_equilibrium::siesta_get_ne_af
real(rprec) function, dimension(:), pointer siesta_get_ne_af(this)
Gets the electron density profile af array.
Definition: siesta_equilibrium.f:1810
vacuum_equilibrium::vacuum_get_b_vec
real(rprec) function, dimension(3) vacuum_get_b_vec(this, x_cart, cyl)
Gets the magnetic field vector at a position.
Definition: vacuum_equilibrium.f:358
equilibrium::equilibrium_get_sxrem_af
real(rprec) function, dimension(:), pointer equilibrium_get_sxrem_af(this, index)
Get the soft x-ray emissivity profile af array.
Definition: equilibrium.f:1833
equilibrium::equilibrium_get_b_vec
real(rprec) function, dimension(3) equilibrium_get_b_vec(this, x_cart, cyl)
Gets the magnetic field vector at a position.
Definition: equilibrium.f:2238
siesta_equilibrium::siesta_get_grid_stop
real(rprec) function siesta_get_grid_stop(this)
Get radial grid size.
Definition: siesta_equilibrium.f:3869
siesta_equilibrium::siesta_get_area_int_element
real(rprec) function siesta_get_area_int_element(this)
Get area integration element.
Definition: siesta_equilibrium.f:3447
siesta_equilibrium::siesta_is_scaler_value
logical function siesta_is_scaler_value(this, id)
Checks if a parameter id is a scaler value.
Definition: siesta_equilibrium.f:3899
siesta_equilibrium::siesta_save_state
subroutine siesta_save_state(this)
Save the internal state of the equilibrium.
Definition: siesta_equilibrium.f:4247
equilibrium::equilibrium_get_magnetic_volume_jphigrid
real(rprec) function, dimension(:,:,:), pointer equilibrium_get_magnetic_volume_jphigrid(this)
Get volume magnetic volume integration jphi grid points.
Definition: equilibrium.f:2542
equilibrium::equilibrium_get_ne
Interface for the equilibrium density profile values.
Definition: equilibrium.f:91
equilibrium::equilibrium_get_gp_sxrem_num_hyper_param
integer function equilibrium_get_gp_sxrem_num_hyper_param(this, index)
Get the number of soft x-ray emissivity gp kernel hyper parameters.
Definition: equilibrium.f:1786
siesta_equilibrium::siesta_def_result
subroutine siesta_def_result(this, result_ncid, maxnsetps_dim_id)
Define NetCDF variables for the result file.
Definition: siesta_equilibrium.f:4506
equilibrium::equilibrium_is_recon_param
logical function equilibrium_is_recon_param(this, id)
Checks if a parameter id is a reconstruction parameter.
Definition: equilibrium.f:3259
equilibrium::equilibrium_get_p_cart
real(rprec) function equilibrium_get_p_cart(this, x_cart)
Gets the plasma pressure at a cartesian position.
Definition: equilibrium.f:2142
equilibrium::equilibrium_get_gp_ti_pp
real(rprec) function equilibrium_get_gp_ti_pp(this, x_cart, y_cart)
Get the ion temperature gp kernel value for the position and position.
Definition: equilibrium.f:1550
siesta_equilibrium::siesta_get_te
Interface for the siesta electron temperature profile values.
Definition: siesta_equilibrium.f:293
equilibrium::equilibrium_write_step_data
subroutine equilibrium_write_step_data(this, result_ncid, current_step)
Write step data to NetCDF result file.
Definition: equilibrium.f:3702
siesta_equilibrium::siesta_is_using_point
logical function siesta_is_using_point(this)
Checks if a the point magnetics are being used.
Definition: siesta_equilibrium.f:4060
equilibrium::equilibrium_siesta_type
integer, parameter equilibrium_siesta_type
Type descriptor for the equilibrium_class subclass siesta_equilibrium.
Definition: equilibrium.f:38
equilibrium::equilibrium_get_gp_sxrem
Interface for the equilibrium guassian process soft x-ray emissivity profile values.
Definition: equilibrium.f:161
siesta_equilibrium::siesta_get_con_surface_krgrid
real(rprec) function, dimension(:,:), pointer siesta_get_con_surface_krgrid(this)
Get the conducting surface integration kr grid points.
Definition: siesta_equilibrium.f:3351
siesta_equilibrium::siesta_get_b_vec
real(rprec) function, dimension(3) siesta_get_b_vec(this, x_cart, cyl)
Gets the magnetic field vector at a position.
Definition: siesta_equilibrium.f:2986
equilibrium
Defines the base class of the type equilibrium_class. Equilibrium is an abstract interface to the equ...
Definition: equilibrium.f:19
equilibrium::equilibrium_construct_vacuum
type(equilibrium_class) function, pointer equilibrium_construct_vacuum(vacuum_object, force_solve)
Construct a equilibrium_class containing a vacuum_equilibrium object.
Definition: equilibrium.f:227
equilibrium::equilibrium_restart
subroutine equilibrium_restart(this, result_ncid, current_step)
Restart equilibrium.
Definition: equilibrium.f:3746
equilibrium::equilibrium_get_gp_ne
Interface for the equilibrium guassian process density profile values.
Definition: equilibrium.f:99
equilibrium::equilibrium_get_ti_af
real(rprec) function, dimension(:), pointer equilibrium_get_ti_af(this)
Get the ion temperature profile af array.
Definition: equilibrium.f:1402
equilibrium::equilibrium_get_param_id
integer function equilibrium_get_param_id(this, param_name)
Get the id for a reconstruction parameter.
Definition: equilibrium.f:532
equilibrium::equilibrium_get_te
Interface for the equilibrium electron temperature profile values.
Definition: equilibrium.f:108
equilibrium::equilibrium_get_gp_ne_pp
real(rprec) function equilibrium_get_gp_ne_pp(this, x_cart, y_cart)
Get the electron density gp kernel value for the position and position.
Definition: equilibrium.f:881
equilibrium::equilibrium_get_gp_sxrem_pp
real(rprec) function equilibrium_get_gp_sxrem_pp(this, x_cart, y_cart, index)
Get the soft x-ray emissivity gp kernel value for the position and position.
Definition: equilibrium.f:1989
equilibrium::equilibrium_get_ti_radial
real(rprec) function equilibrium_get_ti_radial(this, r)
Gets the ion temperature at a radial position.
Definition: equilibrium.f:1647
siesta_equilibrium::siesta_get_magnetic_volume_jphigrid
real(rprec) function, dimension(:,:,:), pointer siesta_get_magnetic_volume_jphigrid(this)
Get volume magnetic volume integration jphi grid points.
Definition: siesta_equilibrium.f:3253
equilibrium::equilibrium_get_gp_te_ij
real(rprec) function equilibrium_get_gp_te_ij(this, i, j)
Get the electron temperature gp kernel value for the two indicies.
Definition: equilibrium.f:1113
siesta_equilibrium::siesta_get_plasma_edge
integer function siesta_get_plasma_edge(this, phi, r, z)
Gets the r and z positions of the outer surface at a toroidal angle.
Definition: siesta_equilibrium.f:3128
equilibrium::equilibrium_converge
logical function equilibrium_converge(this, num_iter, iou, eq_comm, state_flags)
Solves the equilibrium.
Definition: equilibrium.f:3362
siesta_equilibrium::siesta_get_gp_te_num_hyper_param
integer function siesta_get_gp_te_num_hyper_param(this)
Get the number of electron temperature gp kernel hyper parameters.
Definition: siesta_equilibrium.f:2026
siesta_equilibrium::siesta_set_magnetic_cache
Interface for the setting the siesta_magnetic_cache types using siesta_set_magnetic_cache_responce,...
Definition: siesta_equilibrium.f:267
siesta_equilibrium::siesta_converge
logical function siesta_converge(this, num_iter, iou, eq_comm, state_flags)
Solves the siesta equilibrium.
Definition: siesta_equilibrium.f:4107
equilibrium::equilibrium_sync_child
subroutine equilibrium_sync_child(this, index, recon_comm)
Syncronize a child equilibrium state to the parent.
Definition: equilibrium.f:3828
equilibrium::equilibrium_get_p_radial
real(rprec) function equilibrium_get_p_radial(this, r)
Gets the plasma pressure at a radial position.
Definition: equilibrium.f:2189
siesta_equilibrium::siesta_get_grid_size
integer function siesta_get_grid_size(this)
Get radial grid size.
Definition: siesta_equilibrium.f:3811
equilibrium::equilibrium_vacuum_type
integer, parameter equilibrium_vacuum_type
Type descriptor for the equilibrium_class subclass vacuum_equilibrium.
Definition: equilibrium.f:35
equilibrium::equilibrium_set_magnetic_cache
Interface for the equilibrium set magnetic cache.
Definition: equilibrium.f:83
equilibrium::equilibrium_get_con_surface_kphigrid
real(rprec) function, dimension(:,:), pointer equilibrium_get_con_surface_kphigrid(this)
Get the conducting surface integration kphi grid points.
Definition: equilibrium.f:2725
siesta_equilibrium::siesta_get_param_id
integer function siesta_get_param_id(this, param_name)
Get the id for a reconstruction parameter.
Definition: siesta_equilibrium.f:1527
equilibrium::equilibrium_get_te_radial
real(rprec) function equilibrium_get_te_radial(this, r)
Gets the electron temperature at a radial position.
Definition: equilibrium.f:1313
vacuum_equilibrium::vacuum_get_ext_currents
real(rprec) function, dimension(:), pointer vacuum_get_ext_currents(this, scale_currents)
Get external current.
Definition: vacuum_equilibrium.f:457
siesta_equilibrium::siesta_get_ti
Interface for the siesta ion temperature profile values.
Definition: siesta_equilibrium.f:312
vacuum_equilibrium::vacuum_get_int_b_dphi
real(rprec) function vacuum_get_int_b_dphi(this, r, theta)
Gets the loop integrated magnetic field at a position.
Definition: vacuum_equilibrium.f:404
equilibrium::equilibrium_construct
Interface for the construction of equilibrium_class types using equilibrium_construct_vmec,...
Definition: equilibrium.f:74
vacuum_equilibrium::vacuum_destruct
subroutine vacuum_destruct(this)
Deconstruct a vacuum_class object.
Definition: vacuum_equilibrium.f:166
siesta_equilibrium::siesta_get_ti_af
real(rprec) function, dimension(:), pointer siesta_get_ti_af(this)
Gets the ion temperature profile af array.
Definition: siesta_equilibrium.f:2317
siesta_equilibrium::siesta_reset_state
subroutine siesta_reset_state(this)
Reset the internal state of the equilibrium.
Definition: siesta_equilibrium.f:4300
equilibrium::equilibrium_set_magnetic_cache_point
subroutine equilibrium_set_magnetic_cache_point(this, use_axi, state_flags)
Sets the magnetic cache of the equilibrium for points.
Definition: equilibrium.f:475
equilibrium::equilibrium_get_ne_radial
real(rprec) function equilibrium_get_ne_radial(this, r)
Gets the electron density at a radial position.
Definition: equilibrium.f:978
equilibrium::equilibrium_get_ti
Interface for the equilibrium ion temperature profile values.
Definition: equilibrium.f:126
equilibrium::equilibrium_get_gp_ti_ij
real(rprec) function equilibrium_get_gp_ti_ij(this, i, j)
Get the ion temperature gp kernel value for the two indicies.
Definition: equilibrium.f:1448
equilibrium::equilibrium_get_gp_sxrem_ij
real(rprec) function equilibrium_get_gp_sxrem_ij(this, i, j, index)
Get the electron density gp kernel value for the two indicies.
Definition: equilibrium.f:1881
equilibrium::equilibrium_save_state
subroutine equilibrium_save_state(this)
Save the internal state of the equilibrium.
Definition: equilibrium.f:3449
equilibrium::equilibrium_get_magnetic_volume_zgrid
real(rprec) function, dimension(:,:,:), pointer equilibrium_get_magnetic_volume_zgrid(this)
Get volume magnetic volume integration z grid points.
Definition: equilibrium.f:2450
equilibrium::equilibrium_get_con_surface_kzgrid
real(rprec) function, dimension(:,:), pointer equilibrium_get_con_surface_kzgrid(this)
Get the conducting surface integration kz grid points.
Definition: equilibrium.f:2771
equilibrium::equilibrium_set_param
subroutine equilibrium_set_param(this, id, i_index, j_index, value, eq_comm, state_flags)
Sets the value of a reconstruction equilibrium parameter.
Definition: equilibrium.f:353
siesta_equilibrium::siesta_write_init_data
subroutine siesta_write_init_data(this, result_ncid)
Write inital data to NetCDF result file.
Definition: siesta_equilibrium.f:4535
equilibrium::equilibrium_get_gp_ti_pi
real(rprec) function equilibrium_get_gp_ti_pi(this, x_cart, i)
Get the ion temperature gp kernel value for the position and index.
Definition: equilibrium.f:1498
siesta_equilibrium::siesta_write_input
subroutine siesta_write_input(this, current_step)
Write the current valid input.
Definition: siesta_equilibrium.f:4424
siesta_equilibrium::siesta_get_ext_currents
real(rprec) function, dimension(:), pointer siesta_get_ext_currents(this, num_currents, scale_currents)
Get external current.
Definition: siesta_equilibrium.f:3483
equilibrium::equilibrium_get_ext_b_plasma
real(rprec) function, dimension(3) equilibrium_get_ext_b_plasma(this, position, axi_only)
Get external plasma magnetic field.
Definition: equilibrium.f:2932
equilibrium::equilibrium_destruct
subroutine equilibrium_destruct(this)
Deconstruct a equilibrium_class object.
Definition: equilibrium.f:300
siesta_equilibrium::siesta_is_1d_array
logical function siesta_is_1d_array(this, id)
Checks if a parameter id is a 1d array.
Definition: siesta_equilibrium.f:3936
equilibrium::equilibrium_get_te_af
real(rprec) function, dimension(:), pointer equilibrium_get_te_af(this)
Get the electron temperature profile af array.
Definition: equilibrium.f:1067
model_state
Contains parameters defining the bit positions for flags that mark changes in different parts of the ...
Definition: model_state.f:11
siesta_equilibrium::siesta_get_magnetic_volume_jrgrid
real(rprec) function, dimension(:,:,:), pointer siesta_get_magnetic_volume_jrgrid(this)
Get volume magnetic volume integration jr grid points.
Definition: siesta_equilibrium.f:3223
equilibrium::equilibrium_get_ze_radial
real(rprec) function equilibrium_get_ze_radial(this, r)
Gets the effective charge at a radial position.
Definition: equilibrium.f:1740
siesta_equilibrium::siesta_get_ext_b_plasma
real(rprec) function, dimension(3) siesta_get_ext_b_plasma(this, position, axi_only)
Get external plasma magnetic field.
Definition: siesta_equilibrium.f:3521
siesta_equilibrium::siesta_get_te_af
real(rprec) function, dimension(:), pointer siesta_get_te_af(this)
Gets the electron temperature profile af array.
Definition: siesta_equilibrium.f:2061
siesta_equilibrium::siesta_get_sxrem
Interface for the siesta soft x-ray emission profile values.
Definition: siesta_equilibrium.f:329
siesta_equilibrium::siesta_get_p
Interface for the siesta pressure profile values.
Definition: siesta_equilibrium.f:347
equilibrium::equilibrium_get_grid_size
integer function equilibrium_get_grid_size(this)
Get radial grid size.
Definition: equilibrium.f:2979
vacuum_equilibrium::vacuum_get_param_value
real(rprec) function vacuum_get_param_value(this, id, i_index)
Get the value of a reconstruction vacuum parameter.
Definition: vacuum_equilibrium.f:279
equilibrium::equilibrium_is_scaler_value
logical function equilibrium_is_scaler_value(this, id)
Checks if a parameter id is a scaler value.
Definition: equilibrium.f:3115
equilibrium::equilibrium_get_ze_cart
real(rprec) function equilibrium_get_ze_cart(this, x_cart)
Gets the effective charge at a cartesian position.
Definition: equilibrium.f:1693
equilibrium::equilibrium_get_gp_te
Interface for the equilibrium guassian process electron temperature profile values.
Definition: equilibrium.f:117
siesta_equilibrium::siesta_write_step_data
subroutine siesta_write_step_data(this, result_ncid, current_step)
Write step data to NetCDF result file.
Definition: siesta_equilibrium.f:4567
equilibrium::equilibrium_vmec_type
integer, parameter equilibrium_vmec_type
Type descriptor for the equilibrium_class subclass vmec_equilibrium.
Definition: equilibrium.f:32
siesta_equilibrium::siesta_destruct
subroutine siesta_destruct(this)
Deconstruct a siesta_class object.
Definition: siesta_equilibrium.f:496
siesta_equilibrium::siesta_get_param_value
real(rprec) function siesta_get_param_value(this, id, i_index, j_index)
Get the value of a reconstruction siesta parameter.
Definition: siesta_equilibrium.f:1610
vacuum_equilibrium::vacuum_write
subroutine vacuum_write(this, iou)
Write out the equilibrium to an output file.
Definition: vacuum_equilibrium.f:570
equilibrium::equilibrium_get_magnetic_volume_jzgrid
real(rprec) function, dimension(:,:,:), pointer equilibrium_get_magnetic_volume_jzgrid(this)
Get volume magnetic volume integration jz grid points.
Definition: equilibrium.f:2588
model_state::model_state_vmec_flag
integer, parameter model_state_vmec_flag
VMEC Equilibrium changed bit position.
Definition: model_state.f:26
siesta_equilibrium::siesta_is_2d_array
logical function siesta_is_2d_array(this, id)
Checks if a parameter id is a 2d array.
Definition: siesta_equilibrium.f:3977
vacuum_equilibrium::vacuum_get_param_name
character(len=data_name_length) function vacuum_get_param_name(this, id)
Get the name of a reconstruction vacuum parameter.
Definition: vacuum_equilibrium.f:318
equilibrium::equilibrium_get_ti_cart
real(rprec) function equilibrium_get_ti_cart(this, x_cart)
Gets the ion temperature at a cartesian position.
Definition: equilibrium.f:1600
equilibrium::equilibrium_get_ne_cart
real(rprec) function equilibrium_get_ne_cart(this, x_cart)
Gets the electron density at a cartesian position.
Definition: equilibrium.f:931
siesta_equilibrium::siesta_get_gp_ne
Interface for the siesta guassian process density profile values.
Definition: siesta_equilibrium.f:284
equilibrium::equilibrium_get_sxrem_cart
real(rprec) function equilibrium_get_sxrem_cart(this, x_cart, index)
Gets the soft x-ray emissivity at a cartesian position.
Definition: equilibrium.f:2041
equilibrium::equilibrium_get_con_surface_krgrid
real(rprec) function, dimension(:,:), pointer equilibrium_get_con_surface_krgrid(this)
Get the conducting surface integration kr grid points.
Definition: equilibrium.f:2679
equilibrium::equilibrium_sync_state
subroutine equilibrium_sync_state(this, recon_comm)
Syncronize the equilibrium state to children.
Definition: equilibrium.f:3787
siesta_equilibrium::siesta_get_ne
Interface for the siesta density profile values.
Definition: siesta_equilibrium.f:276
equilibrium::equilibrium_get_gp_ti
Interface for the equilibrium guassian process ion temperature profile values.
Definition: equilibrium.f:135
equilibrium::equilibrium_get_grid_start
real(rprec) function equilibrium_get_grid_start(this)
Get start of the radial grid.
Definition: equilibrium.f:3023
equilibrium::equilibrium_get_int_b_dphi
real(rprec) function equilibrium_get_int_b_dphi(this, r, theta)
Gets the loop integrated magnetic field at a position.
Definition: equilibrium.f:2292
equilibrium::equilibrium_get_gp_sxrem_pi
real(rprec) function equilibrium_get_gp_sxrem_pi(this, x_cart, i, index)
Get the soft x-ray emissivity gp kernel value for the position and index.
Definition: equilibrium.f:1935
siesta_equilibrium::siesta_get_gp_sxrem
Interface for the siesta guassian process soft x-ray emissivity profile values.
Definition: siesta_equilibrium.f:338
equilibrium::equilibrium_read_vac_file
subroutine equilibrium_read_vac_file(this, index, eq_comm)
Loads the vacuum magnetic field file.
Definition: equilibrium.f:3417
equilibrium::equilibrium_get_grid_stop
real(rprec) function equilibrium_get_grid_stop(this)
Get stop of the radial grid.
Definition: equilibrium.f:3067
equilibrium::equilibrium_get_param_value
real(rprec) function equilibrium_get_param_value(this, id, i_index, j_index)
Gets the value of a reconstruction equilibrium parameter.
Definition: equilibrium.f:586
equilibrium::equilibrium_get_magnetic_volume_rgrid
real(rprec) function, dimension(:,:,:), pointer equilibrium_get_magnetic_volume_rgrid(this)
Get volume magnetic volume integration radial grid points.
Definition: equilibrium.f:2404
equilibrium::equilibrium_get_ext_currents
real(rprec) function, dimension(:), pointer equilibrium_get_ext_currents(this, num_currents, scale_currents)
Get external current.
Definition: equilibrium.f:2873
equilibrium::equilibrium_set_magnetic_cache_response
subroutine equilibrium_set_magnetic_cache_response(this, response_object, state_flags)
Sets the magnetic cache of the equilibrium for the magnetic responce.
Definition: equilibrium.f:414
equilibrium::equilibrium_class
Base class representing a equilibrium.
Definition: equilibrium.f:50
vacuum_equilibrium::vacuum_get_param_id
integer function vacuum_get_param_id(this, param_name)
Get the id for a reconstruction parameter.
Definition: vacuum_equilibrium.f:237
equilibrium::equilibrium_write_input
subroutine equilibrium_write_input(this, current_step)
Write the current valid input.
Definition: equilibrium.f:3570
siesta_equilibrium::siesta_get_int_b_dphi
real(rprec) function siesta_get_int_b_dphi(this, r, theta)
Gets the loop integrated magnetic field at a position.
Definition: siesta_equilibrium.f:3042
siesta_equilibrium::siesta_get_grid_start
real(rprec) function siesta_get_grid_start(this)
Get start of the radial grid.
Definition: siesta_equilibrium.f:3843
siesta_equilibrium::siesta_get_magnetic_volume_zgrid
real(rprec) function, dimension(:,:,:), pointer siesta_get_magnetic_volume_zgrid(this)
Get volume magnetic volume integration z grid points.
Definition: siesta_equilibrium.f:3193
siesta_equilibrium::siesta_read_vac_file
subroutine siesta_read_vac_file(this, index, eq_comm)
Loads the vacuum magnetic field file.
Definition: siesta_equilibrium.f:4226
equilibrium::equilibrium_get_gp_te_pi
real(rprec) function equilibrium_get_gp_te_pi(this, x_cart, i)
Get the electron temperature gp kernel value for the position and index.
Definition: equilibrium.f:1164
siesta_equilibrium::siesta_get_gp_ti
Interface for the siesta guassian process ion temperature profile values.
Definition: siesta_equilibrium.f:320
siesta_equilibrium::siesta_class
Base class representing a siesta_equilibrium.
Definition: siesta_equilibrium.f:221
equilibrium::equilibrium_is_1d_array
logical function equilibrium_is_1d_array(this, id)
Checks if a parameter id is a 1d array.
Definition: equilibrium.f:3163
equilibrium::equilibrium_no_type
integer, parameter equilibrium_no_type
Type descriptor for the equilibrium_class no subclass.
Definition: equilibrium.f:29
equilibrium::equilibrium_get_gp_ne_pi
real(rprec) function equilibrium_get_gp_ne_pi(this, x_cart, i)
Get the electron density gp kernel value for the position and index.
Definition: equilibrium.f:829
siesta_equilibrium::siesta_write
subroutine siesta_write(this, iou)
Write out the equilibrium to an output file.
Definition: siesta_equilibrium.f:4354
equilibrium::equilibrium_get_gp_ne_num_hyper_param
integer function equilibrium_get_gp_ne_num_hyper_param(this)
Get the number of electron density gp kernel hyper parameters.
Definition: equilibrium.f:690
equilibrium::equilibrium_get_volume_int_element
real(rprec) function equilibrium_get_volume_int_element(this)
Get volume integration element.
Definition: equilibrium.f:2634
siesta_equilibrium::siesta_get_gp_sxrem_num_hyper_param
integer function siesta_get_gp_sxrem_num_hyper_param(this, index)
Get the number of soft x-ray emission gp kernel hyper parameters.
Definition: siesta_equilibrium.f:2534
vacuum_equilibrium::vacuum_is_1d_array
logical function vacuum_is_1d_array(this, id)
Checks if a parameter id is a 1d array.
Definition: vacuum_equilibrium.f:491
equilibrium::equilibrium_get_magnetic_volume_jrgrid
real(rprec) function, dimension(:,:,:), pointer equilibrium_get_magnetic_volume_jrgrid(this)
Get volume magnetic volume integration jr grid points.
Definition: equilibrium.f:2496
equilibrium::equilibrium_reset_state
subroutine equilibrium_reset_state(this)
Reset the internal state of the equilibrium.
Definition: equilibrium.f:3487
equilibrium::equilibrium_is_using_point
logical function equilibrium_is_using_point(this)
Checks if a the point magnetics are being used.
Definition: equilibrium.f:3309
magnetic_response
Defines the base class of the type magnetic_response_class.
Definition: magnetic_response.f:11
equilibrium::equilibrium_get_sxrem_radial
real(rprec) function equilibrium_get_sxrem_radial(this, r, index)
Gets the soft x-ray emissivity at a radial position.
Definition: equilibrium.f:2092
equilibrium::equilibrium_get_p
Interface for the equilibrium pressure profile values.
Definition: equilibrium.f:170
siesta_equilibrium::siesta_get_sxrem_af
real(rprec) function, dimension(:), pointer siesta_get_sxrem_af(this, index)
Gets the soft x-ray emissivity profile af array.
Definition: siesta_equilibrium.f:2572
equilibrium::equilibrium_get_gp_te_num_hyper_param
integer function equilibrium_get_gp_te_num_hyper_param(this)
Get the number of electron temperature gp kernel hyper parameters.
Definition: equilibrium.f:1022
equilibrium::equilibrium_get_te_cart
real(rprec) function equilibrium_get_te_cart(this, x_cart)
Gets the electron temperature at a cartesian position.
Definition: equilibrium.f:1266
vacuum_equilibrium::vacuum_write_input
subroutine vacuum_write_input(this, current_step)
Write the current valid input.
Definition: vacuum_equilibrium.f:601
siesta_equilibrium::siesta_get_gp_ne_num_hyper_param
integer function siesta_get_gp_ne_num_hyper_param(this)
Get the number of electron density gp kernel hyper parameters.
Definition: siesta_equilibrium.f:1775
equilibrium::equilibrium_get_gp_ne_ij
real(rprec) function equilibrium_get_gp_ne_ij(this, i, j)
Get the electron density gp kernel value for the two indicies.
Definition: equilibrium.f:779
equilibrium::equilibrium_write
subroutine equilibrium_write(this, iou)
Write out the equilibrium to an output file.
Definition: equilibrium.f:3527
siesta_equilibrium::siesta_get_gp_te
Interface for the siesta guassian process electron temperature profile values.
Definition: siesta_equilibrium.f:302
equilibrium::equilibrium_get_ze
Interface for the effective charge profile values.
Definition: equilibrium.f:144
siesta_equilibrium::siesta_set_param
subroutine siesta_set_param(this, id, i_index, j_index, value, eq_comm, state_flags)
Sets the value of a reconstruction equilibrium parameter.
Definition: siesta_equilibrium.f:658
vacuum_equilibrium::vacuum_set_param
subroutine vacuum_set_param(this, id, i_index, value)
Sets the value of a reconstruction equilibrium parameter.
Definition: vacuum_equilibrium.f:198
siesta_equilibrium
Defines the base class of the type siesta_class. This module contains all the code necessary to inter...
Definition: siesta_equilibrium.f:74
siesta_equilibrium::siesta_sync_state
subroutine siesta_sync_state(this, recon_comm)
Syncronize the siesta state to children.
Definition: siesta_equilibrium.f:4634
equilibrium::equilibrium_construct_siesta
type(equilibrium_class) function, pointer equilibrium_construct_siesta(siesta_object, force_solve)
Construct a equilibrium_class containing a siesta_equilibrium object.
Definition: equilibrium.f:264
siesta_equilibrium::siesta_get_magnetic_volume_jzgrid
real(rprec) function, dimension(:,:,:), pointer siesta_get_magnetic_volume_jzgrid(this)
Get volume magnetic volume integration jz grid points.
Definition: siesta_equilibrium.f:3283
siesta_equilibrium::siesta_get_con_surface_kzgrid
real(rprec) function, dimension(:,:), pointer siesta_get_con_surface_kzgrid(this)
Get the conducting surface integration kz grid points.
Definition: siesta_equilibrium.f:3411
vacuum_equilibrium::vacuum_class
Base class representing a vacuum_equilibrium.
Definition: vacuum_equilibrium.f:65
siesta_equilibrium::siesta_restart
subroutine siesta_restart(this, result_ncid, current_step)
Restart from a result file.
Definition: siesta_equilibrium.f:4601
equilibrium::equilibrium_get_gp_te_pp
real(rprec) function equilibrium_get_gp_te_pp(this, x_cart, y_cart)
Get the electron temperature gp kernel value for the position and position.
Definition: equilibrium.f:1216
siesta_equilibrium::siesta_get_param_name
character(len=data_name_length) function siesta_get_param_name(this, id)
Get the name of a reconstruction siesta parameter.
Definition: siesta_equilibrium.f:1693
siesta_equilibrium::siesta_sync_child
subroutine siesta_sync_child(this, index, recon_comm)
Syncronize a child siesta state to the parent.
Definition: siesta_equilibrium.f:4691
model_state::model_state_siesta_flag
integer, parameter model_state_siesta_flag
SIESTA Equilibrium changed bit position.
Definition: model_state.f:28
equilibrium::equilibrium_construct_vmec
type(equilibrium_class) function, pointer equilibrium_construct_vmec(vmec_object, force_solve)
Construct a equilibrium_class containing a vmec_equilibrium object.
Definition: equilibrium.f:190
equilibrium::equilibrium_write_init_data
subroutine equilibrium_write_init_data(this, result_ncid)
Write inital data to NetCDF result file.
Definition: equilibrium.f:3660
siesta_equilibrium::siesta_get_volume_int_element
real(rprec) function siesta_get_volume_int_element(this)
Get volume integration element.
Definition: siesta_equilibrium.f:3319
equilibrium::equilibrium_get_plasma_edge
integer function equilibrium_get_plasma_edge(this, phi, r, z)
Gets the r and z positions of the outer surface at a toroidal angle.
Definition: equilibrium.f:2354
vacuum_equilibrium
Defines the base class of the type vacuum_class. This module contains all the code necessary to inter...
Definition: vacuum_equilibrium.f:34
equilibrium::equilibrium_get_sxrem
Interface for the equilibrium soft x-ray emissivity profile values.
Definition: equilibrium.f:152
equilibrium::equilibrium_get_ne_af
real(rprec) function, dimension(:), pointer equilibrium_get_ne_af(this)
Get the electron density profile af array.
Definition: equilibrium.f:735
siesta_equilibrium::siesta_get_gp_ti_num_hyper_param
integer function siesta_get_gp_ti_num_hyper_param(this)
Get the number of ion temperature gp kernel hyper parameters.
Definition: siesta_equilibrium.f:2282
siesta_equilibrium::siesta_is_recon_param
logical function siesta_is_recon_param(this, id)
Checks if a parameter id is a reconstruction parameter.
Definition: siesta_equilibrium.f:4016