132                                               center, scale, coef, long_print)
 
  134    real(kind=
rp), 
intent(in) :: center(3)
 
  135    real(kind=
rp), 
intent(in) :: scale
 
  136    character(len=*), 
intent(in) :: zone_name
 
  137    integer, 
intent(in) :: zone_id
 
  138    type(
coef_t), 
target, 
intent(in) :: coef
 
  139    logical, 
intent(in) :: long_print
 
  140    integer :: n_pts, glb_n_pts, ierr
 
  141    character(len=1000) :: log_buf
 
  144    this%zone_id = zone_id
 
  147    this%zone_name = zone_name
 
  148    if (long_print ) 
then 
  149       this%print_format = 
'(I7,E20.10,E20.10,E20.10,E20.10,A)'     
  151       this%print_format = 
'(I7,E13.5,E13.5,E13.5,E13.5,A)'     
  160    call this%bc%init_base(this%coef) 
 
  161    call this%bc%mark_zone(this%case%msh%labeled_zones(this%zone_id))
 
  162    call this%bc%finalize()
 
  163    n_pts = this%bc%msk(0)
 
  164    call this%n1%init(n_pts)
 
  165    call this%n2%init(n_pts)
 
  166    call this%n3%init(n_pts)
 
  167    call this%r1%init(n_pts)
 
  168    call this%r2%init(n_pts)
 
  169    call this%r3%init(n_pts)
 
  170    call this%force1%init(n_pts)
 
  171    call this%force2%init(n_pts)
 
  172    call this%force3%init(n_pts)
 
  173    call this%force4%init(n_pts)
 
  174    call this%force5%init(n_pts)
 
  175    call this%force6%init(n_pts)
 
  176    call this%s11msk%init(n_pts)
 
  177    call this%s22msk%init(n_pts)
 
  178    call this%s33msk%init(n_pts)
 
  179    call this%s12msk%init(n_pts)
 
  180    call this%s13msk%init(n_pts)
 
  181    call this%s23msk%init(n_pts)
 
  182    call this%pmsk%init(n_pts)
 
  184         this%n1%x, this%n2%x, this%n3%x, n_pts)
 
  186         this%u%size(), n_pts)
 
  188         this%u%size(), n_pts)
 
  190         this%u%size(), n_pts)
 
  192    call mpi_allreduce(n_pts, glb_n_pts, 1, &
 
  195    call neko_log%section(
'Force/torque calculation')
 
  196    write(log_buf, 
'(A,I4,A,A)') 
'Zone ', zone_id, 
'   ', trim(zone_name)
 
  198    write(log_buf, 
'(A,I6, I6)') 
'Global number of GLL points in zone: ', glb_n_pts
 
  200    write(log_buf, 
'(A,E15.7,E15.7,E15.7)') 
'Average of zone''s coordinates: ', &
 
  201                                            glsum(this%r1%x, n_pts)/glb_n_pts, &
 
  202                                            glsum(this%r2%x, n_pts)/glb_n_pts, &
 
  203                                            glsum(this%r3%x, n_pts)/glb_n_pts
 
  205    write(log_buf, 
'(A,E15.7,E15.7,E15.7)') 
'Center for torque calculation: ', center
 
  207    write(log_buf, 
'(A,E15.7)') 
'Scale: ', scale
 
  211    call cadd(this%r1%x,-center(1), n_pts)
 
  212    call cadd(this%r2%x,-center(2), n_pts)
 
  213    call cadd(this%r3%x,-center(3), n_pts)
 
  214    if (neko_bcknd_device .eq. 1 .and. n_pts .gt. 0) 
then 
 
  248    real(kind=
rp), 
intent(in) :: t
 
  249    integer, 
intent(in) :: tstep
 
  250    real(kind=
rp) :: dgtq(12) = 0.0_rp
 
  251    integer :: n_pts, temp_indices(6)
 
  252    type(
field_t), 
pointer :: s11, s22, s33, s12, s13, s23
 
  253    character(len=1000) :: log_buf
 
  255    n_pts = this%bc%msk(0)
 
  265         s13%x, s23%x, this%u, this%v, this%w, this%coef)
 
  268    if (neko_bcknd_device .eq. 0) 
then 
  270            this%u%size(), n_pts)
 
  272            this%u%size(), n_pts)
 
  274            this%u%size(), n_pts)
 
  276            this%u%size(), n_pts)
 
  278            this%u%size(), n_pts)
 
  280            this%u%size(), n_pts)
 
  282            this%u%size(), n_pts)
 
  284                             this%force4%x, this%force5%x, this%force6%x, &
 
  295                             this%case%fluid%mu, &
 
  297       dgtq(1) = 
glsum(this%force1%x, n_pts)
 
  298       dgtq(2) = 
glsum(this%force2%x, n_pts)
 
  299       dgtq(3) = 
glsum(this%force3%x, n_pts)
 
  300       dgtq(4) = 
glsum(this%force4%x, n_pts)
 
  301       dgtq(5) = 
glsum(this%force5%x, n_pts)
 
  302       dgtq(6) = 
glsum(this%force6%x, n_pts)
 
  307       call vcross(this%s11msk%x, this%s22msk%x, this%s33msk%x, &
 
  308                   this%r1%x, this%r2%x, this%r3%x, &
 
  309                   this%force1%x, this%force2%x, this%force3%x, n_pts)
 
  311       dgtq(7) = 
glsum(this%s11msk%x, n_pts)
 
  312       dgtq(8) = 
glsum(this%s22msk%x, n_pts)
 
  313       dgtq(9) = 
glsum(this%s33msk%x, n_pts)
 
  317       call vcross(this%s11msk%x, this%s22msk%x, this%s33msk%x, &
 
  318                    this%r1%x, this%r2%x, this%r3%x, &
 
  319                    this%force4%x, this%force5%x, this%force6%x, n_pts)
 
  320       dgtq(10) = 
glsum(this%s11msk%x, n_pts)
 
  321       dgtq(11) = 
glsum(this%s22msk%x, n_pts)
 
  322       dgtq(12) = 
glsum(this%s33msk%x, n_pts)
 
  324       if (n_pts .gt. 0) 
then 
  326                                      this%bc%msk_d, this%u%size(), n_pts)
 
  328                                      this%bc%msk_d, this%u%size(), n_pts)
 
  330                                      this%bc%msk_d, this%u%size(), n_pts)
 
  332                                      this%bc%msk_d, this%u%size(), n_pts)
 
  334                                      this%bc%msk_d, this%u%size(), n_pts)
 
  336                                      this%bc%msk_d, this%u%size(), n_pts)
 
  338                                      this%bc%msk_d, this%u%size(), n_pts)
 
  355                             this%case%fluid%mu, &
 
  360                             this%r1%x_d, this%r2%x_d, this%r3%x_d, &
 
  361                             this%force1%x_d, this%force2%x_d, &
 
  362                             this%force3%x_d, n_pts)
 
  363          call device_vcross(this%s12msk%x_d,this%s13msk%x_d,this%s23msk%x_d, &
 
  364                      this%r1%x_d, this%r2%x_d, this%r3%x_d, &
 
  365                      this%force4%x_d, this%force5%x_d, this%force6%x_d, n_pts)
 
  380    dgtq = this%scale*dgtq
 
  381    write(log_buf,
'(A, I4, A, A)') 
'Force and torque on zone ', &
 
  382          this%zone_id,
'  ', this%zone_name
 
  384    write(log_buf,
'(A)') &
 
  385          'Time step, time, total force/torque, pressure, viscous, direction' 
  387    write(log_buf, this%print_format) &
 
  388          tstep,t,dgtq(1)+dgtq(4),dgtq(1),dgtq(4),
', forcex' 
  390    write(log_buf, this%print_format) &
 
  391         tstep,t,dgtq(2)+dgtq(5),dgtq(2),dgtq(5),
', forcey' 
  393    write(log_buf, this%print_format) &
 
  394         tstep,t,dgtq(3)+dgtq(6),dgtq(3),dgtq(6),
', forcez' 
  396    write(log_buf, this%print_format) &
 
  397         tstep,t,dgtq(7)+dgtq(10),dgtq(7),dgtq(10),
', torquex' 
  399    write(log_buf, this%print_format) &
 
  400         tstep,t,dgtq(8)+dgtq(11),dgtq(8),dgtq(11),
', torquey' 
  402    write(log_buf, this%print_format) &
 
  403         tstep,t,dgtq(9)+dgtq(12),dgtq(9),dgtq(12),
', torquez'