Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Dynamically change of the blade mass properties within a given simulation in OpoenFAST and ElastoDyn #339

Closed
LaurenceWETI opened this issue Oct 1, 2019 · 15 comments

Comments

@LaurenceWETI
Copy link

Dear OpenFAST developer,

@jjonkman

I’m trying to change blade mass properties dynamically within a given simulation #253. For this reason I tried to modify the blade element mass during the simulation by adding Subroutine AddMass to ElastoDyn.f90 as follow:


END SUBROUTINE Coeff
!----------------------------------------------------------------------------------------------------------------------------------
!**********************************************************************************************************************************
!> This routine is used to add mass to blade elements at a given simulation time
SUBROUTINE AddMass(p)
!..................................................................................................................................

      ! Passed variables

   TYPE(ED_ParameterType),        INTENT(INOUT)    :: p                             !< Parameters of the structural dynamics module


      ! Local variables.


   INTEGER(IntKi)               :: J                                               ! Loops through nodes / elements.
   INTEGER(IntKi)               :: K                                               ! Loops through blades.

   !...............................................................................................................................
   ! Calculate the new blade element mass 
   !...............................................................................................................................

  print*,"Hello Laurence u r in Subroutine AddMass"

   DO K = 1,p%NumBl          ! Loop through the blades

  print*,"Hello Laurence u r in Blades Loop"

      DO J = 1,p%BldNodes    ! Loop through the blade nodes / elements


      ! Calculate the mass of the current element

  print*,"Hello Laurence u r in Nodes Loop"

         p%MassB(K,J) = 1.5*p%MassB(K,J)
         p%BElmntMass(J,K) = p%MassB(K,J)*p%DRNodes(J)                        ! Mass of blade element J


      ENDDO ! J - Blade nodes / elements

   ENDDO    ! K - Blades

END SUBROUTINE AddMass

Subroutine AddMass is called at a given time within the simulation in FAST_Prog.f90 as follow:


!...............................................................................................................................
   ! Time Stepping:
   !...............................................................................................................................         
   
   DO n_t_global = Restart_step, Turbine(1)%p_FAST%n_TMax_m1 
   
      ! call AddMass routine at given simulation time 
	  
	  IF (n_t_global == 100) THEN
	  print*,"Hello Laurence u r at simulation time 100"
			CALL AddMass(p_ED)
	  END IF
	

After that I could compile OpenFAST without errors.
To check if the changing of blade element mass takes place, I performed the TEST 18 (5MW_Land_Dll_WTurb) with the new OpenFAST.exe and I compared the results with the original OpenFAST.exe. The results were identical for both models, and there are no remarkable changes in the results at the given time!
To test whether OpenFAST.f90 is calling Subroutine AddMass correctly or whether the DO Loop in Subroutine AddMass is working properly, some pop-up messages are printed when Subroutines are calling.
cmd

The above screenshot shows clearly that the program didn’t go through the DO Loop. I’m not sure if the ElastoDyn parameters type (p) are available for Subroutine AddMass at the given time or not?
When the ElastoDyn parameters type (p) are defined at initialization and remains stored throughout the entire simulation, how can I assess this parameter within the simulation? And whether it’s possible to change this parameter after initialization?
Any help would be really appreciated!

Best regards

Laurence

@bjonkman
Copy link
Contributor

bjonkman commented Oct 1, 2019

One thing that looks strange to me is that you are using p_ED in FAST_Prog.f90. I would expect it to be Turbine(1)%ED%p. If you added a p_ED type yourself, none of the variables would be initialized, so it either won't go through the loops, or it would give you segmentation faults.

@LaurenceWETI
Copy link
Author

Dear Bonnie,

Thank you very much for your quick response.
I do add a p_ED type by myself in the header of FAST_prog.f90, because I thought that the argument in Fortran Subroutine works in the same way as in Matlab Function, but obviously it seems not.
After I used Turbine(1)%ED%p in FAST_Prog.f90 it did work and it went through the loops and it changed the blade element mass.
Now I’m checking if the Blade Mode Shapes and the Flap- and Edge-stiffness are also changed with the changing of blade element mass.
Could you please keep this issue open? Because I’ll have for sure more questions relating to this topic in the future.

Best regards
Laurence

@LaurenceWETI
Copy link
Author

LaurenceWETI commented Oct 17, 2019

Dear @bjonkman ,

Since the last post I did some modifications to the ElstoDyn.f90 source code, so that not only the blade element mass p%BElmntMass(J,K) is changed within the simulation, but also all the parameters in SUBROUTINE Coeff(p,InputFileData, ErrStat, ErrMsg). To do that I generated a SUBROUTINE UpdateCoeff(p, ErrStat, ErrMsg), which is eventually a copy of Subroutine Coeff() with only one difference – the Subroutine UpdateCoeff() has no Data Type ED_InputFile, because Subroutine UpdateCoeff would be called inside Subroutine ED_CalcOutput mid simulation, whereas the InputFileData is already destroyed after the initialisation. the call of Subroutine UpdateCoeff works as folow:

SUBROUTINE ED_CalcOutput( t, u, p, x, xd, z, OtherState, y, m, ErrStat, ErrMsg )
!..................................................................................................................................
!...............................................................................................................................
   ! set the values in the AllOuts array:
   !...............................................................................................................................
      ! Define the output channel specifying the current simulation time:

   m%AllOuts(  Time) = REAL( t, ReKi )

      ! Change Blade Element Masses:

	IF ( t == 0.50 ) THEN

		p%MassB(K,J) = 1.5*p%MassB(K,J)

		! start with recalculation of the Parameters in subroutine Ceoff()

		CALL UpdateCoeff(p, ErrStat, ErrMsg)
		
	END IF
! This part is done to generate the blade mass element SpnMb as an output to see if the masses changed or not, and it did work  
   DO K = 1,p%NumBl
      DO I = 1, p%NBlGages
		 
				BlElmMass = p%BElmntMass(I,K)
			
		m%AllOuts(  SpnMb(I,K)) = REAL( BlElmMass, R8Ki)

      END DO !I
   END DO !K

But now, because of no Data Type ED_Inputfile in Subroutine UpdateCoeff, I generated new parameters p% to replace the InputFileData% in Subroutine Coeff and to use them after that in Subroutine UpdateCoeff as folw:

SUBROUTINE Coeff(p,InputFileData, ErrStat, ErrMsg)
!..................................................................................................................................
   ELSE                    ! 3-blader
      p%Hubg1Iner = InputFileData%HubIner
      p%HubInerX = InputFileData%HubIner                ! Line 4799 where the 1st error occured 
      p%Hubg2Iner = 0.0
   ENDIF
  
! Calculate the tower shape functions (all derivatives):

      p%TwrFASF(1,J,2) = SHP( p%HNodesNorm(J), p%TwrFlexL, InputFileData%TwFAM1Sh(:), 2, ErrStat, ErrMsg )
      p%TwrFASF(2,J,2) = SHP( p%HNodesNorm(J), p%TwrFlexL, InputFileData%TwFAM2Sh(:), 2, ErrStat, ErrMsg )
      p%TwrFASF(1,J,1) = SHP( p%HNodesNorm(J), p%TwrFlexL, InputFileData%TwFAM1Sh(:), 1, ErrStat, ErrMsg )
      p%TwrFASF(2,J,1) = SHP( p%HNodesNorm(J), p%TwrFlexL, InputFileData%TwFAM2Sh(:), 1, ErrStat, ErrMsg )
      p%TwrFASF(1,J,0) = SHP( p%HNodesNorm(J), p%TwrFlexL, InputFileData%TwFAM1Sh(:), 0, ErrStat, ErrMsg )
      p%TwrFASF(2,J,0) = SHP( p%HNodesNorm(J), p%TwrFlexL, InputFileData%TwFAM2Sh(:), 0, ErrStat, ErrMsg )

	  p%TwFAM1ShX(:) = InputFileData%TwFAM1Sh(:)     ! Line 5195 where the 2nd error occured 
	  p%TwFAM2ShX(:) = InputFileData%TwFAM2Sh(:)

END SUBROUTINE Coeff

!**********************************************************************************************************************************

   SUBROUTINE UpdateCoeff(p, ErrStat, ErrMsg)
!..................................................................................................................................
    ELSE                    ! 3-blader
      p%Hubg1Iner = p%HubInerX
      p%Hubg2Iner = 0.0
   ENDIF

      ! Calculate the tower shape functions (all derivatives):

      p%TwrFASF(1,J,2) = SHP( p%HNodesNorm(J), p%TwrFlexL, p%TwFAM1ShX(:), 2, ErrStat, ErrMsg )

I did that for all InputFileData% in Subroutine Coeff(), and I also allocated the new defined parameter in SetOtherParameters () as folw:

SUBROUTINE SetOtherParameters( p, InputFileData, ErrStat, ErrMsg )
!..................................................................................................................................
  CALL AllocAry(p%TwFAM1ShX,  PolyOrd,			'TwFAM1ShX' , ErrStat, ErrMsg ); IF ( ErrStat /= ErrID_None ) RETURN

But when I start to compile I got this error:

grafik

grafik

I think it’s maybe my defention of the new parameters wrong.
I hope I could explain the problem, and I would be pleased if you could help me to solve this error.

Best regards
Laurence

PS: I took only some sections from the whole code from each Subroutine to illustrate the problem

@bjonkman
Copy link
Contributor

One thing to keep in mind is that if you are changing any registry input files, you might need to build twice (the first time Visual Studio regenerates the *_Types.f90 files, but the code that depends on those types files doesn't necessarily build with the correct version of the *_Types.f90 file until the second build). That may be the cause of the first error you display (unless something else is incorrect).

One issue I see is that your p%TwFAM1SHX array doesn't appear to be the same size as the InputFileData%TwFAM1Sh array. You should allocate them to be the same size if you are going to set them equal to each other. Here's the allocation I see for the existing array:

CALL AllocAry  ( InputFileData%TwFAM1Sh,  PolyOrd-1, 'TwFAM1Sh'  , ErrStat, ErrMsg )

Yours is allocated to be one element larger.

@andrew-platt
Copy link
Collaborator

andrew-platt commented Oct 17, 2019

If you are using the current dev branch, the build process for the _Types.f90 files has changed: they are no longer part of the standard build process. Instead see this: https://raf-openfast.readthedocs.io/en/feature-registry_rework/source/dev/types_files.html

edit: see @rafmudaf comment below for link.

@rafmudaf
Copy link
Collaborator

Actually, this is a better link to the docs: https://openfast.readthedocs.io/en/dev/source/dev/types_files.html.

@LaurenceWETI
Copy link
Author

LaurenceWETI commented Oct 21, 2019

Dear OpenFAST Devloper,
@jjonkman
@bjonkman
I could change the mass of the blade elements within simulation, and I would like to discuss with you the result of the mass addition.
I rised the mass of the blade elements by multiplying the mass density of each blade element with a factor of 1.5. The mass addition proceeded as a part of the Subroutine ED-CalcOutput by calling the Subroutine UpdateCoeff at a given simulation time (in my example at t = 0.5 [s]). I attached the ElastoDyn.f90 file to see all the modifications in my code.

ElastoDyn.txt

After compiling OpenFAST I performed Test 18 (5MW_Land_DLL_WTurb) and I plotted the blade element masses for the List of blade nodes that have strain gages BldGagNd in figure one, this signals (SpnMb) I generated by myself to see the physical-mass addition. In figure 2 I compared the rotor speed of each version of the OpenFAST and of the OpenFAST+AddMass function.

Figure_1.pdf

Figure_2.pdf

Actually, I expiated a continuous reduction of rotor speed after a mass addition of 1.5 for each blade element, but the figures didn’t show that!
Is my expectation right? Did my implementation of the mass addition correct?

Best regards
Laurence

@LaurenceWETI
Copy link
Author

@jjonkman @bjonkman
One more question, I would like to expand the ElastoDyn.dat input file to use the increasing factor
AddMassF and the Time AddMassT for mass addition as InputData, so I don’t have to compile OpenFAST every time I need to change the time or the factor of mass addition:
Plot_1

I already add those input data and theirs corresponding parameters in ElastoDyn_Registry.txt and in ElastoDyn_IO.f90.
I could after that compile OpenFAST, but when I start a simulation I've got this error:
Plot_2

@bjonkman
Copy link
Contributor

Are you sure you are reading the two new inputs from the file correctly? You added two new CALL ReadVar() lines and you put them after the space where it reads PtfmYIner?

@LaurenceWETI
Copy link
Author

Hello,

ist hat posible to use signals from ServoDyn in ElstoDyn? I see that there are signals generated in ElastoDyn to use in ServoDyn like this:

 !...............................................................................................................................
  ! Outputs required for ServoDyn
  !...............................................................................................................................
  
  y%Yaw      = x%QT( DOF_Yaw)
  y%YawRate  = x%QDT(DOF_Yaw)
  y%YawAngle = x%QT( DOF_Yaw) + x%QT(DOF_Y)  !crude approximation for yaw error... (without subtracting it from the wind direction)   
  y%BlPitch  = u%BlPitchCom !OtherState%BlPitch
  y%LSS_Spd  = x%QDT(DOF_GeAz)
  y%HSS_Spd  = ABS(p%GBRatio)*x%QDT(DOF_GeAz)
  y%RotSpeed = x%QDT(DOF_GeAz) + x%QDT(DOF_DrTr)

Now I’m trying to do exact the opposite. I want to generate a signal in ServoDyn to use it in ElasoDyn.
To do that I defined first an output signal
y%BlPitchED(K) = p%BlPitchF(K)
Then I allocated this signal in SrvD_OutputType data structure (y)

   CALL AllocAry( y%BlPitchED, p%NumBl, 'BlPitchED', ErrStat2, ErrMsg2 )
      CALL CheckError( ErrStat2, ErrMsg2 )
      IF (ErrStat >= AbortErrLev) RETURN

Afterthat I defended this output in ServoDyn_Registry.txt
typedef ^ OutputType ReKi BlPitchED {:} - - "Commanded blade pitch angles (ElastoDyn)" radians
But now I don’t know what I’ve to do use this signal in ElstoDyn?
Thank you for your help

@jjonkman
Copy link
Collaborator

Dear @LaurenceWETI,

In the FAST / OpenFAST modularization framework, when modules are coupled together, the inputs and outputs should be consistent, i.e., the inputs to a given module should be derivable from the outputs of other modules. So, if you wish to add an output to ServoDyn that becomes an input to ElastoDyn, you should define the output in ServoDyn, define the input in ElastoDyn, and then modify the OpenFAST glue code to transfer the ServoDyn output as an input to ElastoDyn. For example of this, see how ServoDyn output BlPitchCom is transferred as an input to ElastoDyn.

Best regards,

@LaurenceWETI
Copy link
Author

@jjonkman thank you for your quick response

@daansl1
Copy link

daansl1 commented Jul 14, 2021

Dear all,

@jjonkman
I to have been trying to dynamically adjust parameters in OpenFAST and ElastoDyn, which I have done according to the example @LaurenceWETI so kindly supplied in his first comment in this issue.

I have added the following routine to ElastoDyn.f90 to change GBRatio:

!----------------------------------------------------------------------------------------------------------------------------------
!> This routine changes the gearbox ratio
SUBROUTINE ChangeGBRatio(p)
!..................................................................................................................................

        ! Passed variables

    TYPE(ED_ParameterType),     INTENT(INOUT)   :: p    !< Parameters of the structural dynamics module

        ! Local variables.

    print*, "Changing GBRatio"
    print*, "GBRatio before change:",p%GBRatio
    p%GBRatio = 2*p%GBRatio
    print*, "GBRatio after change:",p%GBRatio

END SUBROUTINE ChangeGBRatio

In order to call the subroutine in the time loop I have added to 'FAST_Prog.f90' the following code:

   !...............................................................................................................................
   ! Time Stepping:
   !...............................................................................................................................         
   
TIME_STEP_LOOP:  DO n_t_global = Restart_step, Turbine(1)%p_FAST%n_TMax_m1 
      
      ! bjj: we have to make sure the n_TMax_m1 and n_ChkptTime are the same for all turbines or have some different logic here
      IF (n_t_global == 6400) THEN 
          print*, "Printing is successfull", n_t_global
          
            CALL ChangeGBRatio(Turbine(1)%ED%p)
      END IF 

Which works perfectly fine when I run OpenFAST from the command line, i.e. when I change the GBRatio, the GenSpeed changes immediately. However, I am using OpenFAST in a Simulink environment and the dynamic change does not occur when I use the S-function.

I have tried calling 'ChangeGBRatio' from the subroutine FAST_Update(iTurb, NumInputs_c, NumOutputs_c, InputAry, OutputAry, ErrStat_c, ErrMsg_c) in 'FAST_library.f90' of the OpenFAST_Simulink-project, similarly to method used for the regular project. The project builds and the S-function is generated without any errors, however the change does not seem to happen during my Simulink run.

REAL(ReKi)                            :: Outputs(NumOutputs_c-1)
   INTEGER(IntKi)                        :: i
   INTEGER(IntKi)                        :: ErrStat2                                ! Error status
   CHARACTER(IntfStrLen-1)               :: ErrMsg2                                 ! Error message  (this needs to be static so that it will print in Matlab's mex library)
                 
   
   IF (n_t_global == 6400) THEN 
          print*, "Printing is successfull", n_t_global
          
            CALL ChangeGBRatio(Turbine(1)%ED%p)
   END IF 
   
   IF ( n_t_global > Turbine(iTurb)%p_FAST%n_TMax_m1 ) THEN !finish 

How can I, similarly to solution for the regular OpenFAST-project, dynamically change the parameter?

If you would be able to help me or point me in the right direction that would be much appreciated. Thank you in advance.

Regards,

Daan

@jjonkman
Copy link
Collaborator

Dear @daansl1 ,

I see that you are using Turbine(1)%ED%p in your CALL to SUBROUTINE ChangeGBRatio(), but turbine index iTurb is set to 0 (being zero-based) in FAST_SFunc.c. Does your code work as expected if you change to Turbine(iTurb)%ED%p?

Best regards,

@daansl1
Copy link

daansl1 commented Jul 15, 2021

Dear @jjonkman ,

Thank you for your swift response. The change you suggested has done the trick, the parameter changes midsimulation. Thank you very much.

Best regards,

Sign up for free to join this conversation on GitHub. Already have an account? Sign in to comment
Projects
None yet
Development

No branches or pull requests

6 participants