V3FIT
utilities.f90
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 !
9 !*******************************************************************************
10  MODULE utilities
11  USE stel_kinds
12  USE island_params, ONLY: ohs=>ohs_i, ns=>ns_i, mpol=>mpol_i, &
13  ntor=>ntor_i, nfp=>nfp_i
14  USE v3_utilities, ONLY: assert, assert_eq
15  USE fourier, ONLY: m0, m1, m2, n0, n1, f_sin, f_cos
16 
17  IMPLICIT NONE
18 
19 !*******************************************************************************
20 ! INTERFACE BLOCKS
21 !*******************************************************************************
22 !-------------------------------------------------------------------------------
25 !-------------------------------------------------------------------------------
26  INTERFACE to_full_mesh
27  MODULE PROCEDURE to_full_mesh2, to_full_mesh3
28  END INTERFACE
29 
30 !-------------------------------------------------------------------------------
33 !-------------------------------------------------------------------------------
34  INTERFACE to_half_mesh
35  MODULE PROCEDURE to_half_mesh, to_half_mesh_p
36  END INTERFACE
37 
38 !-------------------------------------------------------------------------------
41 !-------------------------------------------------------------------------------
42  INTERFACE gradienthalf
43  MODULE PROCEDURE gradienthalf, gradienthalf_p
44  END INTERFACE
45 
46 !-------------------------------------------------------------------------------
49 !-------------------------------------------------------------------------------
50  INTERFACE curl_htof
51  MODULE PROCEDURE curl_htof, curl_htof_targ
52  END INTERFACE
53 
54 ! FIXME: Move to fourier
55 !-------------------------------------------------------------------------------
59 !-------------------------------------------------------------------------------
60  INTERFACE set_bndy_half_to_full
61  MODULE PROCEDURE set_bndy_half_to_full, set_bndy_half_to_full_vec
62  END INTERFACE
63 
64 ! FIXME: Move to fourier
65 !-------------------------------------------------------------------------------
69 !-------------------------------------------------------------------------------
70  INTERFACE set_bndy_half_to_full_ds
71  MODULE PROCEDURE set_bndy_half_to_full_ds, set_bndy_half_to_full_uv_ds
72  END INTERFACE
73 
74 ! FIXME: Move to fourier
75 !-------------------------------------------------------------------------------
78 !-------------------------------------------------------------------------------
79  INTERFACE set_bndy_fouier_m0
80  MODULE PROCEDURE set_bndy_fouier_m0, set_bndy_fouier_m0_vec
81  END INTERFACE
82 
83  CONTAINS
84 
85 !*******************************************************************************
86 ! UTILITY SUBROUTINES
87 !*******************************************************************************
88 !-------------------------------------------------------------------------------
95 !-------------------------------------------------------------------------------
96  SUBROUTINE gradienthalf(gradienth, vecf)
97 
98  IMPLICIT NONE
99 
100 ! Declare Arguments
101  REAL (dp), DIMENSION(:,:,:), ALLOCATABLE, INTENT(inout) :: gradienth
102  REAL (dp), DIMENSION(:,:,:), ALLOCATABLE, INTENT(in) :: vecf
103 
104 ! local variables
105  INTEGER :: nsmin
106  INTEGER :: nsmax
107 
108 ! Start of executable code
109  nsmin = lbound(gradienth,3)
110  nsmax = ubound(gradienth,3)
111 
112  gradienth(:,:,nsmin + 1:nsmax) = ohs*(vecf(:,:,nsmin + 1:nsmax) &
113  - vecf(:,:,nsmin:nsmax - 1))
114  gradienth(:,:,nsmin) = 0
115 
116  END SUBROUTINE
117 
118 !-------------------------------------------------------------------------------
126 !-------------------------------------------------------------------------------
127  SUBROUTINE gradienthalf_p(gradienth, vecf)
128 
129  IMPLICIT NONE
130 
131 ! Declare Arguments
132  REAL (dp), DIMENSION(:,:,:), ALLOCATABLE, INTENT(inout) :: gradienth
133  REAL (dp), DIMENSION(:,:,:), POINTER :: vecf
134 
135 ! local variables
136  INTEGER :: nsmin
137  INTEGER :: nsmax
138 
139 ! Start of executable code
140  nsmin = lbound(gradienth,3)
141  nsmax = ubound(gradienth,3)
142 
143  gradienth(:,:,nsmin + 1:nsmax) = ohs*(vecf(:,:,nsmin + 1:nsmax) &
144  - vecf(:,:,nsmin:nsmax - 1))
145  gradienth(:,:,nsmin) = 0
146 
147  END SUBROUTINE
148 
149 !-------------------------------------------------------------------------------
156 !-------------------------------------------------------------------------------
157  SUBROUTINE gradientfull(gradientf, vech)
158 
159  IMPLICIT NONE
160 
161 ! Declare Arguments
162  REAL (dp), DIMENSION(:,:,:), ALLOCATABLE, INTENT(inout) :: gradientf
163  REAL (dp), DIMENSION(:,:,:), ALLOCATABLE, INTENT(in) :: vech
164 
165 ! local variables
166  INTEGER :: nsmin
167  INTEGER :: nsmax
168 
169 ! Start of executable code
170  nsmin = lbound(vech,3)
171  nsmax = ubound(vech,3)
172 
173  gradientf(:,:,nsmin:nsmax - 1) = ohs*(vech(:,:,nsmin + 1:nsmax) &
174  - vech(:,:,nsmin:nsmax - 1))
175 
176  IF (ubound(gradientf,3) .ge. nsmax) THEN
177  gradientf(:,:,nsmax:) = 0
178  END IF
179 
180  IF (nsmin .eq. 1) THEN
181  gradientf(:,:,1) = 0
182  END IF
183 
184  IF (nsmax .eq. ns) THEN
185  gradientf(:,:,ns) = gradientf(:,:,ns - 1)
186  END IF
187 
188  END SUBROUTINE
189 
190 !-------------------------------------------------------------------------------
197 !-------------------------------------------------------------------------------
198  SUBROUTINE to_half_mesh(vecf, vech)
199 
200  IMPLICIT NONE
201 
202 ! Declare Arguments
203  REAL (dp), DIMENSION(:,:,:), ALLOCATABLE, INTENT(in) :: vecf
204  REAL (dp), DIMENSION(:,:,:), ALLOCATABLE, INTENT(inout) :: vech
205 
206 ! local variables
207  INTEGER :: nsmin
208  INTEGER :: nsmax
209 
210 ! Start of executable code
211  nsmin = lbound(vech,3)
212  nsmax = ubound(vech,3)
213 
214  vech(:,:,nsmin + 1:nsmax) = (vecf(:,:,nsmin:nsmax - 1) &
215  + vecf(:,:,nsmin + 1:nsmax))*0.5
216  vech(:,:,nsmin) = 0
217 
218  END SUBROUTINE
219 
220 !-------------------------------------------------------------------------------
228 !-------------------------------------------------------------------------------
229  SUBROUTINE to_half_mesh_p(vecf, vech)
230 
231  IMPLICIT NONE
232 
233 ! Declare Arguments
234  REAL (dp), DIMENSION(:,:,:), POINTER :: vecf
235  REAL (dp), DIMENSION(:,:,:), ALLOCATABLE, INTENT(inout) :: vech
236 
237 ! local variables
238  INTEGER :: nsmin
239  INTEGER :: nsmax
240 
241 ! Start of executable code
242  nsmin = lbound(vech,3)
243  nsmax = ubound(vech,3)
244  vech(:,:,nsmin + 1:nsmax) = (vecf(:,:,nsmin:nsmax - 1) &
245  + vecf(:,:,nsmin + 1:nsmax))*0.5
246  vech(:,:,nsmin) = 0
247 
248  END SUBROUTINE
249 
250 !-------------------------------------------------------------------------------
259 !-------------------------------------------------------------------------------
260  SUBROUTINE to_full_mesh3(vech, vecf)
261 
262  IMPLICIT NONE
263 
264 ! Declare Arguments
265  REAL (dp), DIMENSION(:,:,:), ALLOCATABLE, INTENT(in) :: vech
266  REAL (dp), DIMENSION(:,:,:), ALLOCATABLE, INTENT(inout) :: vecf
267 
268 ! local variables
269  INTEGER :: nsmin
270  INTEGER :: nsmax
271 
272 ! Start of executable code
273  nsmin = lbound(vecf, 3)
274  CALL assert(nsmin .ge. lbound(vech,3), 'LBOUND WRONG IN to_full_mesh3')
275  nsmax = ubound(vecf, 3)
276  CALL assert(nsmax .le. ubound(vech,3), 'UBOUND WRONG IN to_full_mesh3')
277 
278  vecf(:,:,nsmin:nsmax - 1) = 0.5*(vech(:,:,nsmin + 1:nsmax) &
279  + vech(:,:,nsmin:nsmax - 1))
280  vecf(:,:,nsmax) = 0
281  IF (nsmin .eq. 1) THEN
282  vecf(:,:,1) = vech(:,:,2)
283  END IF
284  IF (nsmax .eq. ns) THEN
285  vecf(:,:,ns) = vech(:,:,ns)
286  END IF
287 
288  END SUBROUTINE
289 
290 !-------------------------------------------------------------------------------
299 !-------------------------------------------------------------------------------
300  SUBROUTINE to_full_mesh2(vech, vecf)
301 
302  IMPLICIT NONE
303 
304 ! Declare Arguments
305  REAL (dp), DIMENSION(:,:), ALLOCATABLE, INTENT(in) :: vech
306  REAL (dp), DIMENSION(:,:), ALLOCATABLE, INTENT(inout) :: vecf
307 
308 ! local variables
309  INTEGER :: nsmin
310  INTEGER :: nsmax
311 
312 ! Start of executable code
313  nsmin = lbound(vecf,2)
314  CALL assert(nsmin .ge. lbound(vech,2), 'LBOUND WRONG IN to_full_mesh2')
315  nsmax = ubound(vecf,2)
316  CALL assert(nsmax .le. ubound(vech,2), 'UBOUND WRONG IN to_full_mesh2')
317 
318  vecf(:,nsmin:nsmax - 1) = 0.5*(vech(:,nsmin:nsmax - 1) &
319  + vech(:,nsmin + 1:nsmax))
320 
321  vecf(:,nsmax) = 0
322  IF (nsmin .EQ. 1) THEN
323  vecf(:,1) = vech(:,2)
324  END IF
325  IF (nsmax .EQ. ns) THEN
326  vecf(:,ns)= vech(:,ns)
327  END IF
328 
329  END SUBROUTINE
330 
331 ! FIXME: Move to a new curl module.
332 !-------------------------------------------------------------------------------
355 !-------------------------------------------------------------------------------
356  SUBROUTINE curl_ftoh(asubsmnf, asubumnf, asubvmnf, &
357  jbsupsmnh, jbsupumnh, jbsupvmnh, parity, &
358  nsmin, nsmax)
359  USE nscalingtools, ONLY: startglobrow
360 
361  IMPLICIT NONE
362 
363 ! Declare Arguments
364  REAL (dp), DIMENSION(:,:,:), POINTER :: asubsmnf
365  REAL (dp), DIMENSION(:,:,:), POINTER :: asubumnf
366  REAL (dp), DIMENSION(:,:,:), POINTER :: asubvmnf
367  REAL (dp), DIMENSION(:,:,:), ALLOCATABLE, INTENT(inout) :: jbsupsmnh
368  REAL (dp), DIMENSION(:,:,:), ALLOCATABLE, INTENT(inout) :: jbsupumnh
369  REAL (dp), DIMENSION(:,:,:), ALLOCATABLE, INTENT(inout) :: jbsupvmnh
370  INTEGER, INTENT(IN) :: parity
371  INTEGER, INTENT(IN) :: nsmin
372  INTEGER, INTENT(IN) :: nsmax
373 
374 ! Local variables
375  INTEGER :: m
376  INTEGER :: n
377  INTEGER :: mp
378  INTEGER :: np
379  INTEGER :: sparity
380  INTEGER :: fouruv
381  INTEGER :: istat
382  INTEGER :: nmin
383  REAL (dp), ALLOCATABLE, DIMENSION(:,:,:) :: asubsmnh
384  REAL (dp), ALLOCATABLE, DIMENSION(:,:,:) :: asubumnh
385  REAL (dp), ALLOCATABLE, DIMENSION(:,:,:) :: asubvmnh
386  REAL (dp), ALLOCATABLE, DIMENSION(:,:,:) :: da_uds
387  REAL (dp), ALLOCATABLE, DIMENSION(:,:,:) :: da_vds
388 
389 ! Start of executable code
390  CALL assert(nsmax .le. ubound(jbsupsmnh,3), 'UBOUND wrong in curl_ftoh')
391 
392  IF (parity .EQ. f_sin) THEN
393  sparity = 1
394  fouruv = f_cos
395  ELSE
396  sparity = -1
397  fouruv = f_sin
398  END IF
399 
400  ALLOCATE(asubsmnh(0:mpol,-ntor:ntor,nsmin:nsmax), &
401  asubumnh(0:mpol,-ntor:ntor,nsmin:nsmax), &
402  asubvmnh(0:mpol,-ntor:ntor,nsmin:nsmax), stat=istat)
403  CALL assert_eq(0, istat, 'Allocation1 failed in curl_ftoh')
404 
405  CALL to_half_mesh(asubsmnf, asubsmnh)
406  CALL to_half_mesh(asubumnf, asubumnh)
407  CALL to_half_mesh(asubvmnf, asubvmnh)
408 
409  CALL set_bndy_fouier_m0(asubsmnh, asubumnh, asubvmnh, parity)
410 
411 ! Compute sqrt(g)B^X contravariant B-field components in Fourier space
412  ALLOCATE(da_uds(0:mpol,-ntor:ntor,nsmin:nsmax), &
413  da_vds(0:mpol,-ntor:ntor,nsmin:nsmax), stat=istat)
414  CALL assert_eq(0, istat, 'Allocation2 failed in curl_ftoh')
415 
416 ! Radial gradients of asubmnu,v on half grid excluding origin and edge. Origin
417 ! and edge will be treated separately.
418  CALL gradienthalf(da_uds, asubumnf)
419  CALL gradienthalf(da_vds, asubvmnf)
420 
421  CALL set_bndy_fouier_m0(da_uds, fouruv)
422  CALL set_bndy_fouier_m0(da_vds, fouruv)
423 
424  nmin = max(2, startglobrow)
425  CALL assert(nmin .ge. lbound(jbsupsmnh,3), 'LBOUND wrong in curl_ftoh')
426 
427 ! (sqrt(g)*B^X)mn
428  DO m = 0, mpol
429  mp = sparity*m
430  DO n = -ntor, ntor
431  np = sparity*n*nfp
432  jbsupsmnh(m,n,nmin:nsmax) = np*asubumnh(m,n,nmin:nsmax) &
433  - mp*asubvmnh(m,n,nmin:nsmax)
434  jbsupumnh(m,n,nmin:nsmax) = np*asubsmnh(m,n,nmin:nsmax) &
435  - da_vds(m,n,nmin:nsmax)
436  jbsupvmnh(m,n,nmin:nsmax) = da_uds(m,n,nmin:nsmax) &
437  - mp*asubsmnh(m,n,nmin:nsmax)
438  END DO
439  END DO
440 
441  DEALLOCATE (asubsmnh, asubumnh, asubvmnh, da_vds, da_uds, stat=istat)
442  CALL assert_eq(0, istat, 'Deallocation failed in curl_ftoh')
443 
444  IF (startglobrow .eq. 1) THEN
445  jbsupsmnh(:,:,1) = 0
446  jbsupumnh(:,:,1) = 0
447  jbsupvmnh(:,:,1) = 0
448  END IF
449 
450  CALL set_bndy_fouier_m0(jbsupsmnh, jbsupumnh, jbsupvmnh, parity)
451 
452  END SUBROUTINE
453 
454 ! FIXME: Move to a new curl module.
455 !-------------------------------------------------------------------------------
480 !-------------------------------------------------------------------------------
481  SUBROUTINE curl_htof(asubsmnh, asubumnh, asubvmnh, &
482  jbsupsmnf, jbsupumnf, jbsupvmnf, parity, &
483  nsmin, nsmax, extend, curtor)
484  USE nscalingtools, ONLY: endglobrow
485  USE stel_constants
486 
487  IMPLICIT NONE
488 
489 ! Declare Arguments
490  REAL (dp), DIMENSION(:,:,:), ALLOCATABLE, INTENT(in) :: asubsmnh
491  REAL (dp), DIMENSION(:,:,:), ALLOCATABLE, INTENT(in) :: asubumnh
492  REAL (dp), DIMENSION(:,:,:), ALLOCATABLE, INTENT(in) :: asubvmnh
493  REAL (dp), DIMENSION(:,:,:), POINTER :: jbsupsmnf
494  REAL (dp), DIMENSION(:,:,:), POINTER :: jbsupumnf
495  REAL (dp), DIMENSION(:,:,:), POINTER :: jbsupvmnf
496  INTEGER, INTENT(in) :: parity
497  INTEGER, INTENT(in) :: nsmin
498  INTEGER, INTENT(in) :: nsmax
499  LOGICAL, INTENT(in) :: extend
500  REAL (dp), INTENT(inout) :: curtor
501 
502 ! Local variables
503  INTEGER :: m
504  INTEGER :: n
505  INTEGER :: mp
506  INTEGER :: np
507  INTEGER :: mo
508  INTEGER :: no
509  INTEGER :: moff
510  INTEGER :: noff
511  INTEGER :: sparity
512  INTEGER :: istat
513  INTEGER :: nmax
514  REAL (dp), DIMENSION(:,:,:), ALLOCATABLE :: asubsmnf
515  REAL (dp), DIMENSION(:,:,:), ALLOCATABLE :: asubumnf
516  REAL (dp), DIMENSION(:,:,:), ALLOCATABLE :: asubvmnf
517  REAL (dp), DIMENSION(:,:,:), ALLOCATABLE :: da_uds
518  REAL (dp), DIMENSION(:,:,:), ALLOCATABLE :: da_vds
519  INTEGER :: fours
520  INTEGER :: fouruv
521 
522 ! Start of executable code
523  CALL assert(nsmin .ge. lbound(jbsupsmnf,3), 'LBOUND wrong in curl_htof')
524 
525  IF (parity .eq. f_sin) THEN
526  sparity = 1
527  fours = f_sin
528  fouruv = f_cos
529  ELSE
530  sparity = -1
531  fours = f_cos
532  fouruv = f_sin
533  END IF
534 
535  moff = lbound(jbsupsmnf,1)
536  noff = lbound(jbsupsmnf,2) + ntor
537 
538  ALLOCATE(asubsmnf(0:mpol,-ntor:ntor,nsmin:nsmax), &
539  asubumnf(0:mpol,-ntor:ntor,nsmin:nsmax), &
540  asubvmnf(0:mpol,-ntor:ntor,nsmin:nsmax), stat=istat)
541  CALL assert_eq(0, istat, 'Allocation1 failed in curl_htof')
542 
543  CALL to_full_mesh(asubsmnh, asubsmnf)
544  CALL to_full_mesh(asubumnh, asubumnf)
545  CALL to_full_mesh(asubvmnh, asubvmnf)
546 
547 ! FIXME: CHECK Boundary Conditions. asubs,subu,subv don't match older version.
548  CALL set_bndy_half_to_full(asubsmnh, asubvmnh, nsmin, parity, &
549  asubsmnf, asubumnf, asubvmnf)
550 
551 ! Store the total enclosed current to check against the orginal VMEC current.
552  IF (parity .eq. f_sin .and. nsmax .eq. ns) THEN
553  curtor = twopi*asubumnf(m0,n0,ns)/mu0
554  END IF
555 
556 ! Compute sqrt(g)B^X == A^X contravariant current components in Fourier space.
557  ALLOCATE(da_uds(0:mpol,-ntor:ntor,nsmin:nsmax), &
558  da_vds(0:mpol,-ntor:ntor,nsmin:nsmax), stat=istat)
559  CALL assert_eq(0, istat, 'Allocation2 failed in curl_htof')
560 
561 ! Radial gradients of asubmnu,v on full grid excluding origin and edge. Origin
562 ! and edge will be treated separately.
563  CALL gradientfull(da_uds, asubumnh)
564  CALL gradientfull(da_vds, asubvmnh)
565 
566  CALL set_bndy_half_to_full_uv_ds(asubumnh, asubvmnh, &
567  da_uds, da_vds, nsmin, fouruv)
568 
569  IF (extend) THEN
570  nmax = min(ns, endglobrow + 1)
571  ELSE
572  nmax = min(ns, endglobrow)
573  END IF
574  CALL assert(nmax .le. ubound(jbsupsmnf,3), 'UBOUND wrong in curl_htof')
575 
576 ! (sqrt(g)*J^X)mn
577  DO m = 0, mpol
578  mp = sparity*m
579  mo = m + moff
580  DO n = -ntor, ntor
581  np = sparity*n*nfp
582  no = n + noff
583  jbsupsmnf(mo,no,nsmin:nmax) = np*asubumnf(m,n,nsmin:nmax) &
584  - mp*asubvmnf(m,n,nsmin:nmax)
585  jbsupumnf(mo,no,nsmin:nmax) = np*asubsmnf(m,n,nsmin:nmax) &
586  - da_vds(m,n,nsmin:nmax)
587  jbsupvmnf(mo,no,nsmin:nmax) = da_uds(m,n,nsmin:nmax) &
588  - mp*asubsmnf(m,n,nsmin:nmax)
589  END DO
590  END DO
591 
592  DEALLOCATE(asubsmnf, asubumnf, asubvmnf, da_uds, da_vds)
593 
594 ! Avoid fp failure in debug mode
595  jbsupsmnf(:,:,nmax+1:ubound(jbsupsmnf,3)) = 0
596  jbsupumnf(:,:,nmax+1:ubound(jbsupsmnf,3)) = 0
597  jbsupvmnf(:,:,nmax+1:ubound(jbsupsmnf,3)) = 0
598 
599  jbsupsmnf(:,:,lbound(jbsupsmnf,3):nsmin-1) = 0
600  jbsupumnf(:,:,lbound(jbsupsmnf,3):nsmin-1) = 0
601  jbsupvmnf(:,:,lbound(jbsupsmnf,3):nsmin-1) = 0
602 
603  CALL set_bndy_fouier_m0(jbsupsmnf, jbsupumnf, jbsupvmnf, parity)
604  IF (nsmin .eq. 1) THEN
605  CALL set_bndy_full_origin(jbsupsmnf, jbsupumnf, jbsupvmnf)
606  END IF
607 
608  END SUBROUTINE
609 
610 ! FIXME: Move to a new curl module.
611 !-------------------------------------------------------------------------------
628 !-------------------------------------------------------------------------------
629  SUBROUTINE curl_htof_targ(asubsmnh, asubumnh, asubvmnh, &
630  jbsupsmnf, jbsupumnf, jbsupvmnf, &
631  parity, nsmin, nsmax, extend, curtor)
632 
633 ! Declare Arguments
634  REAL (dp), ALLOCATABLE, DIMENSION(:,:,:), INTENT(in) :: asubsmnh
635  REAL (dp), ALLOCATABLE, DIMENSION(:,:,:), INTENT(in) :: asubumnh
636  REAL (dp), ALLOCATABLE, DIMENSION(:,:,:), INTENT(in) :: asubvmnh
637  REAL (dp), DIMENSION(:,:,:), ALLOCATABLE, TARGET, INTENT(inout) :: &
638  jbsupsmnf
639  REAL (dp), DIMENSION(:,:,:), ALLOCATABLE, TARGET, INTENT(inout) :: &
640  jbsupumnf
641  REAL (dp), DIMENSION(:,:,:), ALLOCATABLE, TARGET, INTENT(inout) :: &
642  jbsupvmnf
643  INTEGER, INTENT(IN) :: parity
644  INTEGER, INTENT(IN) :: nsmin
645  INTEGER, INTENT(IN) :: nsmax
646  LOGICAL, INTENT(in) :: extend
647  REAL (dp), INTENT(inout) :: curtor
648 
649 ! Local variables
650  REAL (dp), POINTER, DIMENSION(:,:,:) :: jbsp
651  REAL (dp), POINTER, DIMENSION(:,:,:) :: jbup
652  REAL (dp), POINTER, DIMENSION(:,:,:) :: jbvp
653 
654 ! Start of executable code
655  jbsp => jbsupsmnf
656  jbup => jbsupumnf
657  jbvp => jbsupvmnf
658  CALL curl_htof(asubsmnh, asubumnh, asubvmnh, &
659  jbsp, jbup, jbvp, parity, nsmin, nsmax, extend, &
660  curtor)
661 
662  END SUBROUTINE
663 
664 ! FIXME: Move to fourier
665 !-------------------------------------------------------------------------------
677 !-------------------------------------------------------------------------------
678  PURE SUBROUTINE set_bndy_half_to_full(amnh, amnf, nsmin, parity)
679 
680  IMPLICIT NONE
681 
682 ! Declare Arguments
683  REAL (dp), ALLOCATABLE, DIMENSION(:,:,:), INTENT(in) :: amnh
684  REAL (dp), ALLOCATABLE, DIMENSION(:,:,:), INTENT(inout) :: amnf
685  INTEGER, INTENT(IN) :: nsmin
686  INTEGER, INTENT(in) :: parity
687 
688 ! Start of executable code
689  IF (nsmin .eq. 1) THEN
690  amnf(:,:,1) = 0.0
691  amnf(m0:m1,:,1) = amnh(m0:m1,:,2)
692  END IF
693 
694  IF (parity .eq. f_cos) THEN
695  amnf(m0,-ntor:-n1,:) = 0.0
696  ELSE
697  amnf(m0,-ntor:n0,:) = 0.0
698  END IF
699 
700  END SUBROUTINE
701 
702 ! FIXME: Move to fourier
703 !-------------------------------------------------------------------------------
717 !-------------------------------------------------------------------------------
718  PURE SUBROUTINE set_bndy_half_to_full_vec(asubsmnh, asubvmnh, &
719  nsmin, parity, &
720  asubsmnf, asubumnf, asubvmnf)
721 
722  IMPLICIT NONE
723 
724 ! Declare Arguments
725  REAL (dp), ALLOCATABLE, DIMENSION(:,:,:), INTENT(in) :: asubsmnh
726  REAL (dp), ALLOCATABLE, DIMENSION(:,:,:), INTENT(in) :: asubvmnh
727  INTEGER, INTENT(IN) :: nsmin
728  INTEGER, INTENT(in) :: parity
729  REAL (dp), ALLOCATABLE, DIMENSION(:,:,:), INTENT(inout) :: asubsmnf
730  REAL (dp), ALLOCATABLE, DIMENSION(:,:,:), INTENT(inout) :: asubumnf
731  REAL (dp), ALLOCATABLE, DIMENSION(:,:,:), INTENT(inout) :: asubvmnf
732 
733 ! Local Variables
734  INTEGER :: moffh
735  INTEGER :: mofff
736 
737 ! Start of executable code
738  moffh = lbound(asubsmnh, 1)
739  mofff = lbound(asubsmnf, 1)
740 
741  IF (nsmin .eq. 1) THEN
742  asubsmnf(:,:,1) = 0
743  asubsmnf(m1 + mofff,:,1) = asubsmnh(m1 + moffh,:,2)
744 
745  asubumnf(:,:,1) = 0
746 
747  asubvmnf(:,:,1) = 0
748  asubvmnf(m0 + mofff,:,1) = asubvmnh(m0 + moffh,:,2)
749  END IF
750 
751  CALL set_bndy_fouier_m0(asubsmnf, asubumnf, asubvmnf, parity)
752 
753  END SUBROUTINE
754 
755 ! FIXME: Move to fourier
756 !-------------------------------------------------------------------------------
768 !-------------------------------------------------------------------------------
769  PURE SUBROUTINE set_bndy_half_to_full_ds(amnh, amnf_ds, nsmin, parity)
770 
771  IMPLICIT NONE
772 
773 ! Declare Arguments
774  REAL (dp), DIMENSION(:,:,:), INTENT(in) :: amnh
775  REAL (dp), DIMENSION(:,:,:), INTENT(inout) :: amnf_ds
776  INTEGER, INTENT(IN) :: nsmin
777  INTEGER, INTENT(in) :: parity
778 
779 ! Local Variables
780  INTEGER :: moff
781 
782 ! Start of executable code
783  IF (nsmin .eq. 1) THEN
784  moff = lbound(amnh, 1)
785 
786  amnf_ds(:,:,1) = 0.0
787  amnf_ds(m1 + moff,:,1) = 2.0*ohs*amnh(m1 + moff,:,2)
788  END IF
789 
790  CALL set_bndy_fouier_m0(amnf_ds, parity)
791 
792  END SUBROUTINE
793 
794 ! FIXME: Move to fourier
795 !-------------------------------------------------------------------------------
809 !-------------------------------------------------------------------------------
810  PURE SUBROUTINE set_bndy_half_to_full_uv_ds(aumnh, avmnh, &
811  aumnf_ds, avmnf_ds, &
812  nsmin, parity)
813 
814  IMPLICIT NONE
815 
816 ! Declare Arguments
817  REAL (dp), DIMENSION(:,:,:), INTENT(in) :: aumnh
818  REAL (dp), DIMENSION(:,:,:), INTENT(in) :: avmnh
819  REAL (dp), DIMENSION(:,:,:), INTENT(inout) :: aumnf_ds
820  REAL (dp), DIMENSION(:,:,:), INTENT(inout) :: avmnf_ds
821  INTEGER, INTENT(IN) :: nsmin
822  INTEGER, INTENT(in) :: parity
823 
824 ! Local Variables
825  INTEGER :: moff
826 
827 ! Start of executable code
828  IF (nsmin .eq. 1) THEN
829  moff = lbound(aumnh, 1)
830 
831  aumnf_ds(:,:,1) = 0.0
832  aumnf_ds(m0 + moff,:,1) = 2.0*ohs*aumnh(m0 + moff,:,2)
833 
834  avmnf_ds(:,:,1) = 0.0
835  END IF
836 
837  CALL set_bndy_fouier_m0(aumnf_ds, parity)
838  CALL set_bndy_fouier_m0(avmnf_ds, parity)
839 
840  END SUBROUTINE
841 
842 ! FIXME: Move to fourier
843 !-------------------------------------------------------------------------------
853 !-------------------------------------------------------------------------------
854  PURE SUBROUTINE set_bndy_full_origin(amnsf, amnuf, amnvf)
855 
856  IMPLICIT NONE
857 
858 ! Declare Arguments
859  REAL (dp), DIMENSION(:,:,:), INTENT(inout) :: amnsf
860  REAL (dp), DIMENSION(:,:,:), INTENT(inout) :: amnuf
861  REAL (dp), DIMENSION(:,:,:), INTENT(inout) :: amnvf
862 
863 ! Local Variables
864  INTEGER :: moff
865 
866 ! Start of executable code
867  moff = lbound(amnsf, 1)
868 
869  amnsf(m0 + moff,:,1) = 0
870  IF (m2 + moff .lt. ubound(amnsf, 1)) THEN
871  amnsf(m2 + moff:,:,1) = 0
872  END IF
873 
874  amnuf(:,:,1) = 0
875 
876  amnvf(m1 + moff:,:,1) = 0
877 
878  END SUBROUTINE
879 
880 ! FIXME: Move to fourier
881 !-------------------------------------------------------------------------------
889 !-------------------------------------------------------------------------------
890  PURE SUBROUTINE set_bndy_fouier_m0(amn, parity)
891 
892  IMPLICIT NONE
893 
894 ! Declare Arguments
895  REAL (dp), DIMENSION(:,:,:), INTENT(inout) :: amn
896  INTEGER, INTENT(in) :: parity
897 
898 ! Local Variables
899  INTEGER :: moff
900  INTEGER :: noff
901 
902 ! Start of executable code
903  moff = lbound(amn, 1)
904  noff = lbound(amn, 2) + ntor
905 
906  IF (parity .eq. f_cos) THEN
907  amn(m0 + moff,-ntor+noff:-n1+noff,:) = 0
908  ELSE
909  amn(m0 + moff,-ntor+noff:n0+noff,:) = 0
910  END IF
911 
912  END SUBROUTINE
913 
914 ! FIXME: Move to fourier
915 !-------------------------------------------------------------------------------
925 !-------------------------------------------------------------------------------
926  PURE SUBROUTINE set_bndy_fouier_m0_vec(amns, amnu, amnv, parity)
927  USE fourier, ONLY: f_cos, f_sin
928 
929  IMPLICIT NONE
930 
931 ! Declare Arguments
932  REAL (dp), DIMENSION(:,:,:), INTENT(inout) :: amns
933  REAL (dp), DIMENSION(:,:,:), INTENT(inout) :: amnu
934  REAL (dp), DIMENSION(:,:,:), INTENT(inout) :: amnv
935  INTEGER, INTENT(in) :: parity
936 
937 ! Start of executable code
938  IF (parity .eq. f_cos) THEN
939  CALL set_bndy_fouier_m0(amns, f_cos)
940  CALL set_bndy_fouier_m0(amnu, f_sin)
941  CALL set_bndy_fouier_m0(amnv, f_sin)
942  ELSE
943  CALL set_bndy_fouier_m0(amns, f_sin)
944  CALL set_bndy_fouier_m0(amnu, f_cos)
945  CALL set_bndy_fouier_m0(amnv, f_cos)
946  END IF
947 
948  END SUBROUTINE
949 
950  END MODULE
fourier::n1
integer, parameter n1
n = 1 mode.
Definition: fourier.f90:67
v3_utilities::assert_eq
Definition: v3_utilities.f:62
fourier
Module contains subroutines for computing FFT on parallel radial intervals. Converts quantities betwe...
Definition: fourier.f90:13
fourier::m0
integer, parameter m0
m = 0 mode.
Definition: fourier.f90:58
v3_utilities::assert
Definition: v3_utilities.f:55
fourier::f_sin
integer, parameter f_sin
Sine parity.
Definition: fourier.f90:27
island_params
This file contains fix parameters related to the computational grids.
Definition: island_params.f90:10
fourier::n0
integer, parameter n0
n = 0 mode.
Definition: fourier.f90:65
fourier::m2
integer, parameter m2
m = 2 mode.
Definition: fourier.f90:62
fourier::m1
integer, parameter m1
m = 1 mode.
Definition: fourier.f90:60
island_params::ohs_i
real(dp) ohs_i
Radial grid derivative factor. d X/ ds = (x_i+1 - x_i)/(s_i+1 - s_i) (1) Where ohs = 1/(s_i+1 - s_i) ...
Definition: island_params.f90:51
fourier::f_cos
integer, parameter f_cos
Cosine parity.
Definition: fourier.f90:25